1 /* Copyright (c) 2007-2014. The SimGrid Team.
2 * All rights reserved. */
4 /* This program is free software; you can redistribute it and/or modify it
5 * under the terms of the license (GNU LGPL) which comes with this package. */
8 #include "smpi_mpi_dt_private.h"
10 XBT_LOG_NEW_DEFAULT_SUBCATEGORY(smpi_pmpi, smpi,
11 "Logging specific to SMPI (pmpi)");
14 //this function need to be here because of the calls to smpi_bench
15 void TRACE_smpi_set_category(const char *category)
17 //need to end bench otherwise categories for execution tasks are wrong
19 TRACE_internal_smpi_set_category (category);
20 //begin bench after changing process's category
25 /* PMPI User level calls */
27 int PMPI_Init(int *argc, char ***argv)
29 smpi_process_init(argc, argv);
30 smpi_process_mark_as_initialized();
32 int rank = smpi_process_index();
33 TRACE_smpi_init(rank);
34 TRACE_smpi_computing_init(rank);
35 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
36 extra->type = TRACING_INIT;
37 TRACE_smpi_collective_in(rank, -1, __FUNCTION__, extra);
38 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
44 int PMPI_Finalize(void)
48 int rank = smpi_process_index();
49 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
50 extra->type = TRACING_FINALIZE;
51 TRACE_smpi_collective_in(rank, -1, __FUNCTION__, extra);
53 smpi_process_finalize();
55 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
56 TRACE_smpi_finalize(smpi_process_index());
58 smpi_process_destroy();
62 int PMPI_Finalized(int* flag)
64 *flag=smpi_process_finalized();
68 int PMPI_Get_version (int *version,int *subversion){
69 *version = MPI_VERSION;
70 *subversion= MPI_SUBVERSION;
74 int PMPI_Get_library_version (char *version,int *len){
75 int retval = MPI_SUCCESS;
77 snprintf(version,MPI_MAX_LIBRARY_VERSION_STRING,"SMPI Version %d.%d. Copyright The Simgrid Team 2007-2014",SIMGRID_VERSION_MAJOR,
78 SIMGRID_VERSION_MINOR);
79 *len = strlen(version) > MPI_MAX_LIBRARY_VERSION_STRING ? MPI_MAX_LIBRARY_VERSION_STRING : strlen(version);
84 int PMPI_Init_thread(int *argc, char ***argv, int required, int *provided)
86 if (provided != NULL) {
87 *provided = MPI_THREAD_MULTIPLE;
89 return MPI_Init(argc, argv);
92 int PMPI_Query_thread(int *provided)
96 if (provided == NULL) {
99 *provided = MPI_THREAD_MULTIPLE;
100 retval = MPI_SUCCESS;
105 int PMPI_Is_thread_main(int *flag)
110 retval = MPI_ERR_ARG;
112 *flag = smpi_process_index() == 0;
113 retval = MPI_SUCCESS;
118 int PMPI_Abort(MPI_Comm comm, int errorcode)
121 smpi_process_destroy();
122 // FIXME: should kill all processes in comm instead
123 simcall_process_kill(SIMIX_process_self());
127 double PMPI_Wtime(void)
130 if (smpi_process_initialized() && !smpi_process_finalized() && !smpi_process_get_sampling()) {
132 time = SIMIX_get_clock();
135 time = SIMIX_get_clock();
140 extern double sg_maxmin_precision;
141 double PMPI_Wtick(void)
143 return sg_maxmin_precision;
146 int PMPI_Address(void *location, MPI_Aint * address)
151 retval = MPI_ERR_ARG;
153 *address = (MPI_Aint) location;
154 retval = MPI_SUCCESS;
159 int PMPI_Get_address(void *location, MPI_Aint * address)
161 return PMPI_Address(location, address);
164 int PMPI_Type_free(MPI_Datatype * datatype)
169 retval = MPI_ERR_ARG;
171 smpi_datatype_free(datatype);
172 retval = MPI_SUCCESS;
177 int PMPI_Type_size(MPI_Datatype datatype, int *size)
181 if (datatype == MPI_DATATYPE_NULL) {
182 retval = MPI_ERR_TYPE;
183 } else if (size == NULL) {
184 retval = MPI_ERR_ARG;
186 *size = (int) smpi_datatype_size(datatype);
187 retval = MPI_SUCCESS;
192 int PMPI_Type_get_extent(MPI_Datatype datatype, MPI_Aint * lb, MPI_Aint * extent)
196 if (datatype == MPI_DATATYPE_NULL) {
197 retval = MPI_ERR_TYPE;
198 } else if (lb == NULL || extent == NULL) {
199 retval = MPI_ERR_ARG;
201 retval = smpi_datatype_extent(datatype, lb, extent);
206 int PMPI_Type_get_true_extent(MPI_Datatype datatype, MPI_Aint * lb, MPI_Aint * extent)
208 return PMPI_Type_get_extent(datatype, lb, extent);
211 int PMPI_Type_extent(MPI_Datatype datatype, MPI_Aint * extent)
215 if (datatype == MPI_DATATYPE_NULL) {
216 retval = MPI_ERR_TYPE;
217 } else if (extent == NULL) {
218 retval = MPI_ERR_ARG;
220 *extent = smpi_datatype_get_extent(datatype);
221 retval = MPI_SUCCESS;
226 int PMPI_Type_lb(MPI_Datatype datatype, MPI_Aint * disp)
230 if (datatype == MPI_DATATYPE_NULL) {
231 retval = MPI_ERR_TYPE;
232 } else if (disp == NULL) {
233 retval = MPI_ERR_ARG;
235 *disp = smpi_datatype_lb(datatype);
236 retval = MPI_SUCCESS;
241 int PMPI_Type_ub(MPI_Datatype datatype, MPI_Aint * disp)
245 if (datatype == MPI_DATATYPE_NULL) {
246 retval = MPI_ERR_TYPE;
247 } else if (disp == NULL) {
248 retval = MPI_ERR_ARG;
250 *disp = smpi_datatype_ub(datatype);
251 retval = MPI_SUCCESS;
256 int PMPI_Op_create(MPI_User_function * function, int commute, MPI_Op * op)
260 if (function == NULL || op == NULL) {
261 retval = MPI_ERR_ARG;
263 *op = smpi_op_new(function, commute);
264 retval = MPI_SUCCESS;
269 int PMPI_Op_free(MPI_Op * op)
274 retval = MPI_ERR_ARG;
275 } else if (*op == MPI_OP_NULL) {
278 smpi_op_destroy(*op);
280 retval = MPI_SUCCESS;
285 int PMPI_Group_free(MPI_Group * group)
290 retval = MPI_ERR_ARG;
292 smpi_group_destroy(*group);
293 *group = MPI_GROUP_NULL;
294 retval = MPI_SUCCESS;
299 int PMPI_Group_size(MPI_Group group, int *size)
303 if (group == MPI_GROUP_NULL) {
304 retval = MPI_ERR_GROUP;
305 } else if (size == NULL) {
306 retval = MPI_ERR_ARG;
308 *size = smpi_group_size(group);
309 retval = MPI_SUCCESS;
314 int PMPI_Group_rank(MPI_Group group, int *rank)
318 if (group == MPI_GROUP_NULL) {
319 retval = MPI_ERR_GROUP;
320 } else if (rank == NULL) {
321 retval = MPI_ERR_ARG;
323 *rank = smpi_group_rank(group, smpi_process_index());
324 retval = MPI_SUCCESS;
329 int PMPI_Group_translate_ranks(MPI_Group group1, int n, int *ranks1,
330 MPI_Group group2, int *ranks2)
332 int retval, i, index;
333 if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
334 retval = MPI_ERR_GROUP;
336 for (i = 0; i < n; i++) {
337 if(ranks1[i]==MPI_PROC_NULL){
338 ranks2[i]=MPI_PROC_NULL;
340 index = smpi_group_index(group1, ranks1[i]);
341 ranks2[i] = smpi_group_rank(group2, index);
344 retval = MPI_SUCCESS;
349 int PMPI_Group_compare(MPI_Group group1, MPI_Group group2, int *result)
353 if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
354 retval = MPI_ERR_GROUP;
355 } else if (result == NULL) {
356 retval = MPI_ERR_ARG;
358 *result = smpi_group_compare(group1, group2);
359 retval = MPI_SUCCESS;
364 int PMPI_Group_union(MPI_Group group1, MPI_Group group2,
365 MPI_Group * newgroup)
367 int retval, i, proc1, proc2, size, size2;
369 if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
370 retval = MPI_ERR_GROUP;
371 } else if (newgroup == NULL) {
372 retval = MPI_ERR_ARG;
374 size = smpi_group_size(group1);
375 size2 = smpi_group_size(group2);
376 for (i = 0; i < size2; i++) {
377 proc2 = smpi_group_index(group2, i);
378 proc1 = smpi_group_rank(group1, proc2);
379 if (proc1 == MPI_UNDEFINED) {
384 *newgroup = MPI_GROUP_EMPTY;
386 *newgroup = smpi_group_new(size);
387 size2 = smpi_group_size(group1);
388 for (i = 0; i < size2; i++) {
389 proc1 = smpi_group_index(group1, i);
390 smpi_group_set_mapping(*newgroup, proc1, i);
392 for (i = size2; i < size; i++) {
393 proc2 = smpi_group_index(group2, i - size2);
394 smpi_group_set_mapping(*newgroup, proc2, i);
397 retval = MPI_SUCCESS;
402 int PMPI_Group_intersection(MPI_Group group1, MPI_Group group2,
403 MPI_Group * newgroup)
405 int retval, i, proc1, proc2, size;
407 if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
408 retval = MPI_ERR_GROUP;
409 } else if (newgroup == NULL) {
410 retval = MPI_ERR_ARG;
412 size = smpi_group_size(group2);
413 for (i = 0; i < size; i++) {
414 proc2 = smpi_group_index(group2, i);
415 proc1 = smpi_group_rank(group1, proc2);
416 if (proc1 == MPI_UNDEFINED) {
421 *newgroup = MPI_GROUP_EMPTY;
423 *newgroup = smpi_group_new(size);
425 for (i = 0; i < smpi_group_size(group2); i++) {
426 proc2 = smpi_group_index(group2, i);
427 proc1 = smpi_group_rank(group1, proc2);
428 if (proc1 != MPI_UNDEFINED) {
429 smpi_group_set_mapping(*newgroup, proc2, j);
434 retval = MPI_SUCCESS;
439 int PMPI_Group_difference(MPI_Group group1, MPI_Group group2, MPI_Group * newgroup)
441 int retval, i, proc1, proc2, size, size2;
443 if (group1 == MPI_GROUP_NULL || group2 == MPI_GROUP_NULL) {
444 retval = MPI_ERR_GROUP;
445 } else if (newgroup == NULL) {
446 retval = MPI_ERR_ARG;
448 size = size2 = smpi_group_size(group1);
449 for (i = 0; i < size2; i++) {
450 proc1 = smpi_group_index(group1, i);
451 proc2 = smpi_group_rank(group2, proc1);
452 if (proc2 != MPI_UNDEFINED) {
457 *newgroup = MPI_GROUP_EMPTY;
459 *newgroup = smpi_group_new(size);
460 for (i = 0; i < size2; i++) {
461 proc1 = smpi_group_index(group1, i);
462 proc2 = smpi_group_rank(group2, proc1);
463 if (proc2 == MPI_UNDEFINED) {
464 smpi_group_set_mapping(*newgroup, proc1, i);
468 retval = MPI_SUCCESS;
473 int PMPI_Group_incl(MPI_Group group, int n, int *ranks, MPI_Group * newgroup)
475 int retval, i, index;
477 if (group == MPI_GROUP_NULL) {
478 retval = MPI_ERR_GROUP;
479 } else if (newgroup == NULL) {
480 retval = MPI_ERR_ARG;
483 *newgroup = MPI_GROUP_EMPTY;
484 } else if (n == smpi_group_size(group)) {
486 if(group!= smpi_comm_group(MPI_COMM_WORLD)
487 && group != MPI_GROUP_NULL
488 && group != smpi_comm_group(MPI_COMM_SELF)
489 && group != MPI_GROUP_EMPTY)
490 smpi_group_use(group);
492 *newgroup = smpi_group_new(n);
493 for (i = 0; i < n; i++) {
494 index = smpi_group_index(group, ranks[i]);
495 smpi_group_set_mapping(*newgroup, index, i);
498 retval = MPI_SUCCESS;
503 int PMPI_Group_excl(MPI_Group group, int n, int *ranks, MPI_Group * newgroup)
505 int retval, i, j, newsize, oldsize, index;
507 if (group == MPI_GROUP_NULL) {
508 retval = MPI_ERR_GROUP;
509 } else if (newgroup == NULL) {
510 retval = MPI_ERR_ARG;
514 if(group!= smpi_comm_group(MPI_COMM_WORLD)
515 && group != MPI_GROUP_NULL
516 && group != smpi_comm_group(MPI_COMM_SELF)
517 && group != MPI_GROUP_EMPTY)
518 smpi_group_use(group);
519 } else if (n == smpi_group_size(group)) {
520 *newgroup = MPI_GROUP_EMPTY;
522 oldsize=smpi_group_size(group);
523 newsize = oldsize - n;
524 *newgroup = smpi_group_new(newsize);
526 int* to_exclude=xbt_new0(int, smpi_group_size(group));
527 for(i=0; i<oldsize; i++)
530 to_exclude[ranks[i]]=1;
533 for(i=0; i<oldsize; i++){
534 if(to_exclude[i]==0){
535 index = smpi_group_index(group, i);
536 smpi_group_set_mapping(*newgroup, index, j);
541 xbt_free(to_exclude);
543 retval = MPI_SUCCESS;
548 int PMPI_Group_range_incl(MPI_Group group, int n, int ranges[][3],
549 MPI_Group * newgroup)
551 int retval, i, j, rank, size, index;
553 if (group == MPI_GROUP_NULL) {
554 retval = MPI_ERR_GROUP;
555 } else if (newgroup == NULL) {
556 retval = MPI_ERR_ARG;
559 *newgroup = MPI_GROUP_EMPTY;
562 for (i = 0; i < n; i++) {
563 for (rank = ranges[i][0]; /* First */
564 rank >= 0; /* Last */
568 rank += ranges[i][2]; /* Stride */
569 if (ranges[i][0]<ranges[i][1]){
570 if(rank > ranges[i][1])
573 if(rank < ranges[i][1])
579 *newgroup = smpi_group_new(size);
581 for (i = 0; i < n; i++) {
582 for (rank = ranges[i][0]; /* First */
583 rank >= 0; /* Last */
585 index = smpi_group_index(group, rank);
586 smpi_group_set_mapping(*newgroup, index, j);
588 rank += ranges[i][2]; /* Stride */
589 if (ranges[i][0]<ranges[i][1]){
590 if(rank > ranges[i][1])
593 if(rank < ranges[i][1])
599 retval = MPI_SUCCESS;
604 int PMPI_Group_range_excl(MPI_Group group, int n, int ranges[][3],
605 MPI_Group * newgroup)
607 int retval, i, rank, newrank,oldrank, size, index, add;
609 if (group == MPI_GROUP_NULL) {
610 retval = MPI_ERR_GROUP;
611 } else if (newgroup == NULL) {
612 retval = MPI_ERR_ARG;
616 if(group!= smpi_comm_group(MPI_COMM_WORLD)
617 && group != MPI_GROUP_NULL
618 && group != smpi_comm_group(MPI_COMM_SELF)
619 && group != MPI_GROUP_EMPTY)
620 smpi_group_use(group);
622 size = smpi_group_size(group);
623 for (i = 0; i < n; i++) {
624 for (rank = ranges[i][0]; /* First */
625 rank >= 0; /* Last */
629 rank += ranges[i][2]; /* Stride */
630 if (ranges[i][0]<ranges[i][1]){
631 if(rank > ranges[i][1])
634 if(rank < ranges[i][1])
640 *newgroup = MPI_GROUP_EMPTY;
642 *newgroup = smpi_group_new(size);
645 while (newrank < size) {
647 for (i = 0; i < n; i++) {
648 for (rank = ranges[i][0];rank >= 0;){
654 rank += ranges[i][2]; /* Stride */
656 if (ranges[i][0]<ranges[i][1]){
657 if(rank > ranges[i][1])
660 if(rank < ranges[i][1])
666 index = smpi_group_index(group, oldrank);
667 smpi_group_set_mapping(*newgroup, index, newrank);
675 retval = MPI_SUCCESS;
680 int PMPI_Comm_rank(MPI_Comm comm, int *rank)
683 if (comm == MPI_COMM_NULL) {
684 retval = MPI_ERR_COMM;
685 } else if (rank == NULL) {
686 retval = MPI_ERR_ARG;
688 *rank = smpi_comm_rank(comm);
689 retval = MPI_SUCCESS;
694 int PMPI_Comm_size(MPI_Comm comm, int *size)
697 if (comm == MPI_COMM_NULL) {
698 retval = MPI_ERR_COMM;
699 } else if (size == NULL) {
700 retval = MPI_ERR_ARG;
702 *size = smpi_comm_size(comm);
703 retval = MPI_SUCCESS;
708 int PMPI_Comm_get_name (MPI_Comm comm, char* name, int* len)
712 if (comm == MPI_COMM_NULL) {
713 retval = MPI_ERR_COMM;
714 } else if (name == NULL || len == NULL) {
715 retval = MPI_ERR_ARG;
717 smpi_comm_get_name(comm, name, len);
718 retval = MPI_SUCCESS;
723 int PMPI_Comm_group(MPI_Comm comm, MPI_Group * group)
727 if (comm == MPI_COMM_NULL) {
728 retval = MPI_ERR_COMM;
729 } else if (group == NULL) {
730 retval = MPI_ERR_ARG;
732 *group = smpi_comm_group(comm);
733 if(*group!= smpi_comm_group(MPI_COMM_WORLD)
734 && *group != MPI_GROUP_NULL
735 && *group != smpi_comm_group(MPI_COMM_SELF)
736 && *group != MPI_GROUP_EMPTY)
737 smpi_group_use(*group);
738 retval = MPI_SUCCESS;
743 int PMPI_Comm_compare(MPI_Comm comm1, MPI_Comm comm2, int *result)
747 if (comm1 == MPI_COMM_NULL || comm2 == MPI_COMM_NULL) {
748 retval = MPI_ERR_COMM;
749 } else if (result == NULL) {
750 retval = MPI_ERR_ARG;
752 if (comm1 == comm2) { /* Same communicators means same groups */
756 smpi_group_compare(smpi_comm_group(comm1),
757 smpi_comm_group(comm2));
758 if (*result == MPI_IDENT) {
759 *result = MPI_CONGRUENT;
762 retval = MPI_SUCCESS;
767 int PMPI_Comm_dup(MPI_Comm comm, MPI_Comm * newcomm)
771 if (comm == MPI_COMM_NULL) {
772 retval = MPI_ERR_COMM;
773 } else if (newcomm == NULL) {
774 retval = MPI_ERR_ARG;
776 *newcomm = smpi_comm_new(smpi_comm_group(comm), smpi_comm_topo(comm));
777 retval = MPI_SUCCESS;
782 int PMPI_Comm_create(MPI_Comm comm, MPI_Group group, MPI_Comm * newcomm)
786 if (comm == MPI_COMM_NULL) {
787 retval = MPI_ERR_COMM;
788 } else if (group == MPI_GROUP_NULL) {
789 retval = MPI_ERR_GROUP;
790 } else if (newcomm == NULL) {
791 retval = MPI_ERR_ARG;
792 } else if(smpi_group_rank(group,smpi_process_index())==MPI_UNDEFINED){
793 *newcomm= MPI_COMM_NULL;
794 retval = MPI_SUCCESS;
797 *newcomm = smpi_comm_new(group, NULL);
798 retval = MPI_SUCCESS;
803 int PMPI_Comm_free(MPI_Comm * comm)
808 retval = MPI_ERR_ARG;
809 } else if (*comm == MPI_COMM_NULL) {
810 retval = MPI_ERR_COMM;
812 smpi_comm_destroy(*comm);
813 *comm = MPI_COMM_NULL;
814 retval = MPI_SUCCESS;
819 int PMPI_Comm_disconnect(MPI_Comm * comm)
821 /* TODO: wait until all communication in comm are done */
825 retval = MPI_ERR_ARG;
826 } else if (*comm == MPI_COMM_NULL) {
827 retval = MPI_ERR_COMM;
829 smpi_comm_destroy(*comm);
830 *comm = MPI_COMM_NULL;
831 retval = MPI_SUCCESS;
836 int PMPI_Comm_split(MPI_Comm comm, int color, int key, MPI_Comm* comm_out)
841 if (comm_out == NULL) {
842 retval = MPI_ERR_ARG;
843 } else if (comm == MPI_COMM_NULL) {
844 retval = MPI_ERR_COMM;
846 *comm_out = smpi_comm_split(comm, color, key);
847 retval = MPI_SUCCESS;
854 int PMPI_Send_init(void *buf, int count, MPI_Datatype datatype, int dst,
855 int tag, MPI_Comm comm, MPI_Request * request)
860 if (request == NULL) {
861 retval = MPI_ERR_ARG;
862 } else if (comm == MPI_COMM_NULL) {
863 retval = MPI_ERR_COMM;
864 } else if (dst == MPI_PROC_NULL) {
865 retval = MPI_SUCCESS;
867 *request = smpi_mpi_send_init(buf, count, datatype, dst, tag, comm);
868 retval = MPI_SUCCESS;
871 if (retval != MPI_SUCCESS && request)
872 *request = MPI_REQUEST_NULL;
876 int PMPI_Recv_init(void *buf, int count, MPI_Datatype datatype, int src,
877 int tag, MPI_Comm comm, MPI_Request * request)
882 if (request == NULL) {
883 retval = MPI_ERR_ARG;
884 } else if (comm == MPI_COMM_NULL) {
885 retval = MPI_ERR_COMM;
886 } else if (src == MPI_PROC_NULL) {
887 retval = MPI_SUCCESS;
889 *request = smpi_mpi_recv_init(buf, count, datatype, src, tag, comm);
890 retval = MPI_SUCCESS;
893 if (retval != MPI_SUCCESS && request)
894 *request = MPI_REQUEST_NULL;
898 int PMPI_Ssend_init(void* buf, int count, MPI_Datatype datatype,
899 int dst, int tag, MPI_Comm comm, MPI_Request* request)
904 if (request == NULL) {
905 retval = MPI_ERR_ARG;
906 } else if (comm == MPI_COMM_NULL) {
907 retval = MPI_ERR_COMM;
908 } else if (dst == MPI_PROC_NULL) {
909 retval = MPI_SUCCESS;
911 *request = smpi_mpi_ssend_init(buf, count, datatype, dst, tag, comm);
912 retval = MPI_SUCCESS;
915 if (retval != MPI_SUCCESS && request)
916 *request = MPI_REQUEST_NULL;
920 int PMPI_Start(MPI_Request * request)
925 if (request == NULL || *request == MPI_REQUEST_NULL) {
926 retval = MPI_ERR_REQUEST;
928 smpi_mpi_start(*request);
929 retval = MPI_SUCCESS;
935 int PMPI_Startall(int count, MPI_Request * requests)
940 if (requests == NULL) {
941 retval = MPI_ERR_ARG;
943 retval = MPI_SUCCESS;
944 for (i = 0 ; i < count ; i++) {
945 if(requests[i] == MPI_REQUEST_NULL) {
946 retval = MPI_ERR_REQUEST;
949 if(retval != MPI_ERR_REQUEST) {
950 smpi_mpi_startall(count, requests);
957 int PMPI_Request_free(MPI_Request * request)
962 if (*request == MPI_REQUEST_NULL) {
963 retval = MPI_ERR_ARG;
965 smpi_mpi_request_free(request);
966 retval = MPI_SUCCESS;
972 int PMPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src,
973 int tag, MPI_Comm comm, MPI_Request * request)
979 if (request == NULL) {
980 retval = MPI_ERR_ARG;
981 } else if (comm == MPI_COMM_NULL) {
982 retval = MPI_ERR_COMM;
983 } else if (src == MPI_PROC_NULL) {
984 *request = MPI_REQUEST_NULL;
985 retval = MPI_SUCCESS;
986 } else if (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0)){
987 retval = MPI_ERR_COMM;
988 } else if (count < 0) {
989 retval = MPI_ERR_COUNT;
990 } else if (buf==NULL && count > 0) {
991 retval = MPI_ERR_COUNT;
992 } else if (datatype == MPI_DATATYPE_NULL){
993 retval = MPI_ERR_TYPE;
994 } else if(tag<0 && tag != MPI_ANY_TAG){
995 retval = MPI_ERR_TAG;
999 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1000 int src_traced = smpi_group_index(smpi_comm_group(comm), src);
1002 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1003 extra->type = TRACING_IRECV;
1004 extra->send_size = count;
1005 extra->src = src_traced;
1007 extra->datatype1 = encode_datatype(datatype);
1008 TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__, extra);
1011 *request = smpi_mpi_irecv(buf, count, datatype, src, tag, comm);
1012 retval = MPI_SUCCESS;
1015 TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
1016 (*request)->recv = 1;
1021 if (retval != MPI_SUCCESS && request)
1022 *request = MPI_REQUEST_NULL;
1027 int PMPI_Isend(void *buf, int count, MPI_Datatype datatype, int dst,
1028 int tag, MPI_Comm comm, MPI_Request * request)
1033 if (request == NULL) {
1034 retval = MPI_ERR_ARG;
1035 } else if (comm == MPI_COMM_NULL) {
1036 retval = MPI_ERR_COMM;
1037 } else if (dst == MPI_PROC_NULL) {
1038 *request = MPI_REQUEST_NULL;
1039 retval = MPI_SUCCESS;
1040 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1041 retval = MPI_ERR_RANK;
1042 } else if (count < 0) {
1043 retval = MPI_ERR_COUNT;
1044 } else if (buf==NULL && count > 0) {
1045 retval = MPI_ERR_COUNT;
1046 } else if (datatype == MPI_DATATYPE_NULL){
1047 retval = MPI_ERR_TYPE;
1048 } else if(tag<0 && tag != MPI_ANY_TAG){
1049 retval = MPI_ERR_TAG;
1053 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1054 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1056 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1057 extra->type = TRACING_ISEND;
1058 extra->send_size = count;
1060 extra->dst = dst_traced;
1061 extra->datatype1 = encode_datatype(datatype);
1062 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, extra);
1063 TRACE_smpi_send(rank, rank, dst_traced, count*smpi_datatype_size(datatype));
1066 *request = smpi_mpi_isend(buf, count, datatype, dst, tag, comm);
1067 retval = MPI_SUCCESS;
1070 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1071 (*request)->send = 1;
1076 if (retval != MPI_SUCCESS && request)
1077 *request = MPI_REQUEST_NULL;
1081 int PMPI_Issend(void* buf, int count, MPI_Datatype datatype,
1082 int dst, int tag, MPI_Comm comm, MPI_Request* request)
1087 if (request == NULL) {
1088 retval = MPI_ERR_ARG;
1089 } else if (comm == MPI_COMM_NULL) {
1090 retval = MPI_ERR_COMM;
1091 } else if (dst == MPI_PROC_NULL) {
1092 *request = MPI_REQUEST_NULL;
1093 retval = MPI_SUCCESS;
1094 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1095 retval = MPI_ERR_RANK;
1096 } else if (count < 0) {
1097 retval = MPI_ERR_COUNT;
1098 } else if (buf==NULL && count > 0) {
1099 retval = MPI_ERR_COUNT;
1100 } else if (datatype == MPI_DATATYPE_NULL){
1101 retval = MPI_ERR_TYPE;
1102 } else if(tag<0 && tag != MPI_ANY_TAG){
1103 retval = MPI_ERR_TAG;
1107 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1108 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1109 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1110 extra->type = TRACING_ISSEND;
1111 extra->send_size = count;
1113 extra->dst = dst_traced;
1114 extra->datatype1 = encode_datatype(datatype);
1115 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, extra);
1116 TRACE_smpi_send(rank, rank, dst_traced, count*smpi_datatype_size(datatype));
1119 *request = smpi_mpi_issend(buf, count, datatype, dst, tag, comm);
1120 retval = MPI_SUCCESS;
1123 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1124 (*request)->send = 1;
1129 if (retval != MPI_SUCCESS && request)
1130 *request = MPI_REQUEST_NULL;
1134 int PMPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag,
1135 MPI_Comm comm, MPI_Status * status)
1140 if (comm == MPI_COMM_NULL) {
1141 retval = MPI_ERR_COMM;
1142 } else if (src == MPI_PROC_NULL) {
1143 smpi_empty_status(status);
1144 status->MPI_SOURCE = MPI_PROC_NULL;
1145 retval = MPI_SUCCESS;
1146 } else if (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0)){
1147 retval = MPI_ERR_RANK;
1148 } else if (count < 0) {
1149 retval = MPI_ERR_COUNT;
1150 } else if (buf==NULL && count > 0) {
1151 retval = MPI_ERR_COUNT;
1152 } else if (datatype == MPI_DATATYPE_NULL){
1153 retval = MPI_ERR_TYPE;
1154 } else if(tag<0 && tag != MPI_ANY_TAG){
1155 retval = MPI_ERR_TAG;
1158 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1159 int src_traced = smpi_group_index(smpi_comm_group(comm), src);
1160 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1161 extra->type = TRACING_RECV;
1162 extra->send_size = count;
1163 extra->src = src_traced;
1165 extra->datatype1 = encode_datatype(datatype);
1166 TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__, extra);
1169 smpi_mpi_recv(buf, count, datatype, src, tag, comm, status);
1170 retval = MPI_SUCCESS;
1173 //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1174 if(status!=MPI_STATUS_IGNORE){
1175 src_traced = smpi_group_index(smpi_comm_group(comm), status->MPI_SOURCE);
1176 TRACE_smpi_recv(rank, src_traced, rank);
1178 TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
1186 int PMPI_Send(void *buf, int count, MPI_Datatype datatype, int dst, int tag,
1193 if (comm == MPI_COMM_NULL) {
1194 retval = MPI_ERR_COMM;
1195 } else if (dst == MPI_PROC_NULL) {
1196 retval = MPI_SUCCESS;
1197 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1198 retval = MPI_ERR_RANK;
1199 } else if (count < 0) {
1200 retval = MPI_ERR_COUNT;
1201 } else if (buf==NULL && count > 0) {
1202 retval = MPI_ERR_COUNT;
1203 } else if (datatype == MPI_DATATYPE_NULL){
1204 retval = MPI_ERR_TYPE;
1205 } else if(tag<0 && tag != MPI_ANY_TAG){
1206 retval = MPI_ERR_TAG;
1210 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1211 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1212 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1213 extra->type = TRACING_SEND;
1214 extra->send_size = count;
1216 extra->dst = dst_traced;
1217 extra->datatype1 = encode_datatype(datatype);
1218 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, extra);
1219 TRACE_smpi_send(rank, rank, dst_traced,count*smpi_datatype_size(datatype));
1222 smpi_mpi_send(buf, count, datatype, dst, tag, comm);
1223 retval = MPI_SUCCESS;
1226 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1236 int PMPI_Ssend(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm) {
1241 if (comm == MPI_COMM_NULL) {
1242 retval = MPI_ERR_COMM;
1243 } else if (dst == MPI_PROC_NULL) {
1244 retval = MPI_SUCCESS;
1245 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1246 retval = MPI_ERR_RANK;
1247 } else if (count < 0) {
1248 retval = MPI_ERR_COUNT;
1249 } else if (buf==NULL && count > 0) {
1250 retval = MPI_ERR_COUNT;
1251 } else if (datatype == MPI_DATATYPE_NULL){
1252 retval = MPI_ERR_TYPE;
1253 } else if(tag<0 && tag != MPI_ANY_TAG){
1254 retval = MPI_ERR_TAG;
1258 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1259 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1260 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1261 extra->type = TRACING_SSEND;
1262 extra->send_size = count;
1264 extra->dst = dst_traced;
1265 extra->datatype1 = encode_datatype(datatype);
1266 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, extra); TRACE_smpi_send(rank, rank, dst_traced,count*smpi_datatype_size(datatype));
1269 smpi_mpi_ssend(buf, count, datatype, dst, tag, comm);
1270 retval = MPI_SUCCESS;
1273 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1281 int PMPI_Sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1282 int dst, int sendtag, void *recvbuf, int recvcount,
1283 MPI_Datatype recvtype, int src, int recvtag,
1284 MPI_Comm comm, MPI_Status * status)
1290 if (comm == MPI_COMM_NULL) {
1291 retval = MPI_ERR_COMM;
1292 } else if (sendtype == MPI_DATATYPE_NULL
1293 || recvtype == MPI_DATATYPE_NULL) {
1294 retval = MPI_ERR_TYPE;
1295 } else if (src == MPI_PROC_NULL || dst == MPI_PROC_NULL) {
1296 smpi_empty_status(status);
1297 status->MPI_SOURCE = MPI_PROC_NULL;
1298 retval = MPI_SUCCESS;
1299 }else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0 ||
1300 (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0))){
1301 retval = MPI_ERR_RANK;
1302 } else if (sendcount < 0 || recvcount<0) {
1303 retval = MPI_ERR_COUNT;
1304 } else if ((sendbuf==NULL && sendcount > 0)||(recvbuf==NULL && recvcount>0)) {
1305 retval = MPI_ERR_COUNT;
1306 } else if((sendtag<0 && sendtag != MPI_ANY_TAG)||(recvtag<0 && recvtag != MPI_ANY_TAG)){
1307 retval = MPI_ERR_TAG;
1311 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1312 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1313 int src_traced = smpi_group_index(smpi_comm_group(comm), src);
1314 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1315 extra->type = TRACING_SENDRECV;
1316 extra->send_size = sendcount;
1317 extra->recv_size = recvcount;
1318 extra->src = src_traced;
1319 extra->dst = dst_traced;
1320 extra->datatype1 = encode_datatype(sendtype);
1321 extra->datatype2 = encode_datatype(recvtype);
1323 TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__, extra);
1324 TRACE_smpi_send(rank, rank, dst_traced,sendcount*smpi_datatype_size(sendtype));
1328 smpi_mpi_sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf,
1329 recvcount, recvtype, src, recvtag, comm, status);
1330 retval = MPI_SUCCESS;
1333 TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1334 TRACE_smpi_recv(rank, src_traced, rank);
1343 int PMPI_Sendrecv_replace(void *buf, int count, MPI_Datatype datatype,
1344 int dst, int sendtag, int src, int recvtag,
1345 MPI_Comm comm, MPI_Status * status)
1347 //TODO: suboptimal implementation
1350 if (datatype == MPI_DATATYPE_NULL) {
1351 retval = MPI_ERR_TYPE;
1352 } else if (count < 0) {
1353 retval = MPI_ERR_COUNT;
1355 int size = smpi_datatype_get_extent(datatype) * count;
1356 recvbuf = xbt_new0(char, size);
1358 MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf, count,
1359 datatype, src, recvtag, comm, status);
1360 if(retval==MPI_SUCCESS){
1361 smpi_datatype_copy(recvbuf, count, datatype, buf, count, datatype);
1369 int PMPI_Test(MPI_Request * request, int *flag, MPI_Status * status)
1373 if (request == NULL || flag == NULL) {
1374 retval = MPI_ERR_ARG;
1375 } else if (*request == MPI_REQUEST_NULL) {
1377 smpi_empty_status(status);
1378 retval = MPI_ERR_REQUEST;
1381 int rank = request && (*request)->comm != MPI_COMM_NULL
1382 ? smpi_process_index()
1385 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1386 extra->type = TRACING_TEST;
1387 TRACE_smpi_testing_in(rank, extra);
1389 *flag = smpi_mpi_test(request, status);
1391 TRACE_smpi_testing_out(rank);
1393 retval = MPI_SUCCESS;
1399 int PMPI_Testany(int count, MPI_Request requests[], int *index, int *flag,
1400 MPI_Status * status)
1405 if (index == NULL || flag == NULL) {
1406 retval = MPI_ERR_ARG;
1408 *flag = smpi_mpi_testany(count, requests, index, status);
1409 retval = MPI_SUCCESS;
1415 int PMPI_Testall(int count, MPI_Request* requests, int* flag, MPI_Status* statuses)
1421 retval = MPI_ERR_ARG;
1423 *flag = smpi_mpi_testall(count, requests, statuses);
1424 retval = MPI_SUCCESS;
1430 int PMPI_Probe(int source, int tag, MPI_Comm comm, MPI_Status* status) {
1434 if (status == NULL) {
1435 retval = MPI_ERR_ARG;
1436 } else if (comm == MPI_COMM_NULL) {
1437 retval = MPI_ERR_COMM;
1438 } else if (source == MPI_PROC_NULL) {
1439 smpi_empty_status(status);
1440 status->MPI_SOURCE = MPI_PROC_NULL;
1441 retval = MPI_SUCCESS;
1443 smpi_mpi_probe(source, tag, comm, status);
1444 retval = MPI_SUCCESS;
1451 int PMPI_Iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status) {
1456 retval = MPI_ERR_ARG;
1457 } else if (status == NULL) {
1458 retval = MPI_ERR_ARG;
1459 } else if (comm == MPI_COMM_NULL) {
1460 retval = MPI_ERR_COMM;
1461 } else if (source == MPI_PROC_NULL) {
1463 smpi_empty_status(status);
1464 status->MPI_SOURCE = MPI_PROC_NULL;
1465 retval = MPI_SUCCESS;
1467 smpi_mpi_iprobe(source, tag, comm, flag, status);
1468 retval = MPI_SUCCESS;
1474 int PMPI_Wait(MPI_Request * request, MPI_Status * status)
1480 smpi_empty_status(status);
1482 if (request == NULL) {
1483 retval = MPI_ERR_ARG;
1484 } else if (*request == MPI_REQUEST_NULL) {
1485 retval = MPI_ERR_REQUEST;
1489 int rank = request && (*request)->comm != MPI_COMM_NULL
1490 ? smpi_process_index()
1493 int src_traced = (*request)->src;
1494 int dst_traced = (*request)->dst;
1495 MPI_Comm comm = (*request)->comm;
1496 int is_wait_for_receive = (*request)->recv;
1497 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1498 extra->type = TRACING_WAIT;
1499 TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__, extra);
1502 smpi_mpi_wait(request, status);
1503 retval = MPI_SUCCESS;
1506 //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1507 TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1508 if (is_wait_for_receive) {
1509 if(src_traced==MPI_ANY_SOURCE)
1510 src_traced = (status!=MPI_STATUS_IGNORE) ?
1511 smpi_group_rank(smpi_comm_group(comm), status->MPI_SOURCE) :
1513 TRACE_smpi_recv(rank, src_traced, dst_traced);
1523 int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * status)
1530 //save requests information for tracing
1532 int *srcs = xbt_new0(int, count);
1533 int *dsts = xbt_new0(int, count);
1534 int *recvs = xbt_new0(int, count);
1535 MPI_Comm *comms = xbt_new0(MPI_Comm, count);
1537 for (i = 0; i < count; i++) {
1538 MPI_Request req = requests[i]; //already received requests are no longer valid
1542 recvs[i] = req->recv;
1543 comms[i] = req->comm;
1546 int rank_traced = smpi_process_index();
1547 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1548 extra->type = TRACING_WAITANY;
1549 extra->send_size=count;
1550 TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__,extra);
1553 *index = smpi_mpi_waitany(count, requests, status);
1555 if(*index!=MPI_UNDEFINED){
1556 int src_traced = srcs[*index];
1557 //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1558 int dst_traced = dsts[*index];
1559 int is_wait_for_receive = recvs[*index];
1560 if (is_wait_for_receive) {
1561 if(srcs[*index]==MPI_ANY_SOURCE)
1562 src_traced = (status!=MPI_STATUSES_IGNORE) ?
1563 smpi_group_rank(smpi_comm_group(comms[*index]), status->MPI_SOURCE) :
1565 TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1567 TRACE_smpi_ptp_out(rank_traced, src_traced, dst_traced, __FUNCTION__);
1579 int PMPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
1584 //save information from requests
1586 int *srcs = xbt_new0(int, count);
1587 int *dsts = xbt_new0(int, count);
1588 int *recvs = xbt_new0(int, count);
1589 int *valid = xbt_new0(int, count);
1590 MPI_Comm *comms = xbt_new0(MPI_Comm, count);
1592 //int valid_count = 0;
1593 for (i = 0; i < count; i++) {
1594 MPI_Request req = requests[i];
1595 if(req!=MPI_REQUEST_NULL){
1598 recvs[i] = req->recv;
1599 comms[i] = req->comm;
1605 int rank_traced = smpi_process_index();
1606 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1607 extra->type = TRACING_WAITALL;
1608 extra->send_size=count;
1609 TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__,extra);
1611 int retval = smpi_mpi_waitall(count, requests, status);
1613 for (i = 0; i < count; i++) {
1615 //int src_traced = srcs[*index];
1616 //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1617 int src_traced = srcs[i];
1618 int dst_traced = dsts[i];
1619 int is_wait_for_receive = recvs[i];
1620 if (is_wait_for_receive) {
1621 if(src_traced==MPI_ANY_SOURCE)
1622 src_traced = (status!=MPI_STATUSES_IGNORE) ?
1623 smpi_group_rank(smpi_comm_group(comms[i]), status[i].MPI_SOURCE) :
1625 TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1629 TRACE_smpi_ptp_out(rank_traced, -1, -1, __FUNCTION__);
1641 int PMPI_Waitsome(int incount, MPI_Request requests[], int *outcount,
1642 int *indices, MPI_Status status[])
1647 if (outcount == NULL) {
1648 retval = MPI_ERR_ARG;
1650 *outcount = smpi_mpi_waitsome(incount, requests, indices, status);
1651 retval = MPI_SUCCESS;
1657 int PMPI_Testsome(int incount, MPI_Request requests[], int* outcount,
1658 int* indices, MPI_Status status[])
1663 if (outcount == NULL) {
1664 retval = MPI_ERR_ARG;
1666 *outcount = smpi_mpi_testsome(incount, requests, indices, status);
1667 retval = MPI_SUCCESS;
1674 int PMPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm)
1680 if (comm == MPI_COMM_NULL) {
1681 retval = MPI_ERR_COMM;
1684 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1685 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1687 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1688 extra->type = TRACING_BCAST;
1689 extra->send_size = count;
1690 extra->root = root_traced;
1691 extra->datatype1 = encode_datatype(datatype);
1692 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__, extra);
1695 mpi_coll_bcast_fun(buf, count, datatype, root, comm);
1696 retval = MPI_SUCCESS;
1698 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1706 int PMPI_Barrier(MPI_Comm comm)
1712 if (comm == MPI_COMM_NULL) {
1713 retval = MPI_ERR_COMM;
1716 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1717 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1718 extra->type = TRACING_BARRIER;
1719 TRACE_smpi_collective_in(rank, -1, __FUNCTION__, extra);
1721 mpi_coll_barrier_fun(comm);
1722 retval = MPI_SUCCESS;
1724 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1732 int PMPI_Gather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1733 void *recvbuf, int recvcount, MPI_Datatype recvtype,
1734 int root, MPI_Comm comm)
1740 if (comm == MPI_COMM_NULL) {
1741 retval = MPI_ERR_COMM;
1742 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1743 ((smpi_comm_rank(comm) == root) && (recvtype == MPI_DATATYPE_NULL))){
1744 retval = MPI_ERR_TYPE;
1745 } else if ((( sendbuf != MPI_IN_PLACE) && (sendcount <0)) ||
1746 ((smpi_comm_rank(comm) == root) && (recvcount <0))){
1747 retval = MPI_ERR_COUNT;
1750 char* sendtmpbuf = (char*) sendbuf;
1751 int sendtmpcount = sendcount;
1752 MPI_Datatype sendtmptype = sendtype;
1753 if( (smpi_comm_rank(comm) == root) && (sendbuf == MPI_IN_PLACE )) {
1755 sendtmptype=recvtype;
1758 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1759 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1760 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1761 extra->type = TRACING_GATHER;
1762 extra->send_size = sendtmpcount;
1763 extra->recv_size = recvcount;
1764 extra->root = root_traced;
1765 extra->datatype1 = encode_datatype(sendtmptype);
1766 extra->datatype2 = encode_datatype(recvtype);
1768 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__, extra);
1770 mpi_coll_gather_fun(sendtmpbuf, sendtmpcount, sendtmptype, recvbuf, recvcount,
1771 recvtype, root, comm);
1774 retval = MPI_SUCCESS;
1776 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1784 int PMPI_Gatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1785 void *recvbuf, int *recvcounts, int *displs,
1786 MPI_Datatype recvtype, int root, MPI_Comm comm)
1792 if (comm == MPI_COMM_NULL) {
1793 retval = MPI_ERR_COMM;
1794 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1795 ((smpi_comm_rank(comm) == root) && (recvtype == MPI_DATATYPE_NULL))){
1796 retval = MPI_ERR_TYPE;
1797 } else if (( sendbuf != MPI_IN_PLACE) && (sendcount <0)){
1798 retval = MPI_ERR_COUNT;
1799 } else if (recvcounts == NULL || displs == NULL) {
1800 retval = MPI_ERR_ARG;
1802 char* sendtmpbuf = (char*) sendbuf;
1803 int sendtmpcount = sendcount;
1804 MPI_Datatype sendtmptype = sendtype;
1805 if( (smpi_comm_rank(comm) == root) && (sendbuf == MPI_IN_PLACE )) {
1807 sendtmptype=recvtype;
1811 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1812 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1814 int size = smpi_comm_size(comm);
1815 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1816 extra->type = TRACING_GATHERV;
1817 extra->send_size = sendtmpcount;
1818 extra->recvcounts= xbt_malloc(size*sizeof(int));
1819 for(i=0; i< size; i++)//copy data to avoid bad free
1820 extra->recvcounts[i] = recvcounts[i];
1821 extra->num_processes = size;
1822 extra->root = root_traced;
1823 extra->datatype1 = encode_datatype(sendtmptype);
1824 extra->datatype2 = encode_datatype(recvtype);
1826 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__,extra);
1828 smpi_mpi_gatherv(sendtmpbuf, sendtmpcount, sendtmptype, recvbuf, recvcounts,
1829 displs, recvtype, root, comm);
1830 retval = MPI_SUCCESS;
1832 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1840 int PMPI_Allgather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1841 void *recvbuf, int recvcount, MPI_Datatype recvtype,
1848 if (comm == MPI_COMM_NULL) {
1849 retval = MPI_ERR_COMM;
1850 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1851 (recvtype == MPI_DATATYPE_NULL)){
1852 retval = MPI_ERR_TYPE;
1853 } else if ((( sendbuf != MPI_IN_PLACE) && (sendcount <0)) ||
1855 retval = MPI_ERR_COUNT;
1857 if(sendbuf == MPI_IN_PLACE) {
1858 sendbuf=((char*)recvbuf)+smpi_datatype_get_extent(recvtype)*recvcount*smpi_comm_rank(comm);
1859 sendcount=recvcount;
1863 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1864 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1865 extra->type = TRACING_ALLGATHER;
1866 extra->send_size = sendcount;
1867 extra->recv_size = recvcount;
1868 extra->datatype1 = encode_datatype(sendtype);
1869 extra->datatype2 = encode_datatype(recvtype);
1871 TRACE_smpi_collective_in(rank, -1, __FUNCTION__, extra);
1873 mpi_coll_allgather_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1875 retval = MPI_SUCCESS;
1878 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1885 int PMPI_Allgatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1886 void *recvbuf, int *recvcounts, int *displs,
1887 MPI_Datatype recvtype, MPI_Comm comm)
1893 if (comm == MPI_COMM_NULL) {
1894 retval = MPI_ERR_COMM;
1895 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1896 (recvtype == MPI_DATATYPE_NULL)){
1897 retval = MPI_ERR_TYPE;
1898 } else if (( sendbuf != MPI_IN_PLACE) && (sendcount <0)){
1899 retval = MPI_ERR_COUNT;
1900 } else if (recvcounts == NULL || displs == NULL) {
1901 retval = MPI_ERR_ARG;
1904 if(sendbuf == MPI_IN_PLACE) {
1905 sendbuf=((char*)recvbuf)+smpi_datatype_get_extent(recvtype)*displs[smpi_comm_rank(comm)];
1906 sendcount=recvcounts[smpi_comm_rank(comm)];
1910 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1912 int size = smpi_comm_size(comm);
1913 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1914 extra->type = TRACING_ALLGATHERV;
1915 extra->send_size = sendcount;
1916 extra->recvcounts= xbt_malloc(size*sizeof(int));
1917 for(i=0; i< size; i++)//copy data to avoid bad free
1918 extra->recvcounts[i] = recvcounts[i];
1919 extra->num_processes = size;
1920 extra->datatype1 = encode_datatype(sendtype);
1921 extra->datatype2 = encode_datatype(recvtype);
1923 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
1925 mpi_coll_allgatherv_fun(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
1926 displs, recvtype, comm);
1927 retval = MPI_SUCCESS;
1929 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1937 int PMPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1938 void *recvbuf, int recvcount, MPI_Datatype recvtype,
1939 int root, MPI_Comm comm)
1945 if (comm == MPI_COMM_NULL) {
1946 retval = MPI_ERR_COMM;
1947 } else if (((smpi_comm_rank(comm)==root) && (sendtype == MPI_DATATYPE_NULL))
1948 || ((recvbuf !=MPI_IN_PLACE) && (recvtype == MPI_DATATYPE_NULL))) {
1949 retval = MPI_ERR_TYPE;
1952 if (recvbuf == MPI_IN_PLACE) {
1954 recvcount=sendcount;
1957 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1958 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1959 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1960 extra->type = TRACING_SCATTER;
1961 extra->send_size = sendcount;
1962 extra->recv_size= recvcount;
1963 extra->root = root_traced;
1964 extra->datatype1 = encode_datatype(sendtype);
1965 extra->datatype2 = encode_datatype(recvtype);
1967 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__,extra);
1969 mpi_coll_scatter_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1970 recvtype, root, comm);
1971 retval = MPI_SUCCESS;
1973 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1981 int PMPI_Scatterv(void *sendbuf, int *sendcounts, int *displs,
1982 MPI_Datatype sendtype, void *recvbuf, int recvcount,
1983 MPI_Datatype recvtype, int root, MPI_Comm comm)
1989 if (comm == MPI_COMM_NULL) {
1990 retval = MPI_ERR_COMM;
1991 } else if (sendcounts == NULL || displs == NULL) {
1992 retval = MPI_ERR_ARG;
1993 } else if (((smpi_comm_rank(comm)==root) && (sendtype == MPI_DATATYPE_NULL))
1994 || ((recvbuf !=MPI_IN_PLACE) && (recvtype == MPI_DATATYPE_NULL))) {
1995 retval = MPI_ERR_TYPE;
1997 if (recvbuf == MPI_IN_PLACE) {
1999 recvcount=sendcounts[smpi_comm_rank(comm)];
2002 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2003 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
2005 int size = smpi_comm_size(comm);
2006 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2007 extra->type = TRACING_SCATTERV;
2008 extra->recv_size = recvcount;
2009 extra->sendcounts= xbt_malloc(size*sizeof(int));
2010 for(i=0; i< size; i++)//copy data to avoid bad free
2011 extra->sendcounts[i] = sendcounts[i];
2012 extra->num_processes = size;
2013 extra->root = root_traced;
2014 extra->datatype1 = encode_datatype(sendtype);
2015 extra->datatype2 = encode_datatype(recvtype);
2017 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__,extra);
2020 smpi_mpi_scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf,
2021 recvcount, recvtype, root, comm);
2022 retval = MPI_SUCCESS;
2024 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
2032 int PMPI_Reduce(void *sendbuf, void *recvbuf, int count,
2033 MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
2039 if (comm == MPI_COMM_NULL) {
2040 retval = MPI_ERR_COMM;
2041 } else if (datatype == MPI_DATATYPE_NULL || op == MPI_OP_NULL) {
2042 retval = MPI_ERR_ARG;
2045 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2046 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
2047 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2048 extra->type = TRACING_REDUCE;
2049 extra->send_size = count;
2050 extra->datatype1 = encode_datatype(datatype);
2051 extra->root = root_traced;
2053 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__,extra);
2055 mpi_coll_reduce_fun(sendbuf, recvbuf, count, datatype, op, root, comm);
2057 retval = MPI_SUCCESS;
2059 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
2067 int PMPI_Reduce_local(void *inbuf, void *inoutbuf, int count,
2068 MPI_Datatype datatype, MPI_Op op){
2072 if (datatype == MPI_DATATYPE_NULL || op == MPI_OP_NULL) {
2073 retval = MPI_ERR_ARG;
2075 smpi_op_apply(op, inbuf, inoutbuf, &count, &datatype);
2082 int PMPI_Allreduce(void *sendbuf, void *recvbuf, int count,
2083 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2089 if (comm == MPI_COMM_NULL) {
2090 retval = MPI_ERR_COMM;
2091 } else if (datatype == MPI_DATATYPE_NULL) {
2092 retval = MPI_ERR_TYPE;
2093 } else if (op == MPI_OP_NULL) {
2094 retval = MPI_ERR_OP;
2097 char* sendtmpbuf = (char*) sendbuf;
2098 if( sendbuf == MPI_IN_PLACE ) {
2099 sendtmpbuf = (char *)xbt_malloc(count*smpi_datatype_get_extent(datatype));
2100 smpi_datatype_copy(recvbuf, count, datatype,sendtmpbuf, count, datatype);
2103 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2104 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2105 extra->type = TRACING_ALLREDUCE;
2106 extra->send_size = count;
2107 extra->datatype1 = encode_datatype(datatype);
2109 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2111 mpi_coll_allreduce_fun(sendtmpbuf, recvbuf, count, datatype, op, comm);
2113 if( sendbuf == MPI_IN_PLACE ) {
2114 xbt_free(sendtmpbuf);
2117 retval = MPI_SUCCESS;
2119 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2127 int PMPI_Scan(void *sendbuf, void *recvbuf, int count,
2128 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2134 if (comm == MPI_COMM_NULL) {
2135 retval = MPI_ERR_COMM;
2136 } else if (datatype == MPI_DATATYPE_NULL) {
2137 retval = MPI_ERR_TYPE;
2138 } else if (op == MPI_OP_NULL) {
2139 retval = MPI_ERR_OP;
2142 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2143 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2144 extra->type = TRACING_SCAN;
2145 extra->send_size = count;
2146 extra->datatype1 = encode_datatype(datatype);
2148 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2150 smpi_mpi_scan(sendbuf, recvbuf, count, datatype, op, comm);
2151 retval = MPI_SUCCESS;
2153 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2161 int PMPI_Exscan(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype,
2162 MPI_Op op, MPI_Comm comm){
2167 if (comm == MPI_COMM_NULL) {
2168 retval = MPI_ERR_COMM;
2169 } else if (datatype == MPI_DATATYPE_NULL) {
2170 retval = MPI_ERR_TYPE;
2171 } else if (op == MPI_OP_NULL) {
2172 retval = MPI_ERR_OP;
2175 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2176 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2177 extra->type = TRACING_EXSCAN;
2178 extra->send_size = count;
2179 extra->datatype1 = encode_datatype(datatype);
2181 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2183 smpi_mpi_exscan(sendbuf, recvbuf, count, datatype, op, comm);
2184 retval = MPI_SUCCESS;
2186 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2194 int PMPI_Reduce_scatter(void *sendbuf, void *recvbuf, int *recvcounts,
2195 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2200 if (comm == MPI_COMM_NULL) {
2201 retval = MPI_ERR_COMM;
2202 } else if (datatype == MPI_DATATYPE_NULL) {
2203 retval = MPI_ERR_TYPE;
2204 } else if (op == MPI_OP_NULL) {
2205 retval = MPI_ERR_OP;
2206 } else if (recvcounts == NULL) {
2207 retval = MPI_ERR_ARG;
2210 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2212 int size = smpi_comm_size(comm);
2213 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2214 extra->type = TRACING_REDUCE_SCATTER;
2215 extra->send_size = 0;
2216 extra->recvcounts= xbt_malloc(size*sizeof(int));
2217 for(i=0; i< size; i++)//copy data to avoid bad free
2218 extra->recvcounts[i] = recvcounts[i];
2219 extra->num_processes = size;
2220 extra->datatype1 = encode_datatype(datatype);
2222 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2224 void* sendtmpbuf=sendbuf;
2225 if(sendbuf==MPI_IN_PLACE){
2229 mpi_coll_reduce_scatter_fun(sendtmpbuf, recvbuf, recvcounts,
2230 datatype, op, comm);
2231 retval = MPI_SUCCESS;
2233 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2241 int PMPI_Reduce_scatter_block(void *sendbuf, void *recvbuf, int recvcount,
2242 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2247 if (comm == MPI_COMM_NULL) {
2248 retval = MPI_ERR_COMM;
2249 } else if (datatype == MPI_DATATYPE_NULL) {
2250 retval = MPI_ERR_TYPE;
2251 } else if (op == MPI_OP_NULL) {
2252 retval = MPI_ERR_OP;
2253 } else if (recvcount < 0) {
2254 retval = MPI_ERR_ARG;
2256 int count=smpi_comm_size(comm);
2259 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2260 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2261 extra->type = TRACING_REDUCE_SCATTER;
2262 extra->send_size = 0;
2263 extra->recvcounts= xbt_malloc(count*sizeof(int));
2264 for(i=0; i< count; i++)//copy data to avoid bad free
2265 extra->recvcounts[i] = recvcount;
2266 extra->num_processes = count;
2267 extra->datatype1 = encode_datatype(datatype);
2269 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2271 int* recvcounts=(int*)xbt_malloc(count);
2272 for (i=0; i<count;i++)recvcounts[i]=recvcount;
2273 mpi_coll_reduce_scatter_fun(sendbuf, recvbuf, recvcounts,
2274 datatype, op, comm);
2275 xbt_free(recvcounts);
2276 retval = MPI_SUCCESS;
2278 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2286 int PMPI_Alltoall(void *sendbuf, int sendcount, MPI_Datatype sendtype,
2287 void *recvbuf, int recvcount, MPI_Datatype recvtype,
2294 if (comm == MPI_COMM_NULL) {
2295 retval = MPI_ERR_COMM;
2296 } else if (sendtype == MPI_DATATYPE_NULL
2297 || recvtype == MPI_DATATYPE_NULL) {
2298 retval = MPI_ERR_TYPE;
2301 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2302 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2303 extra->type = TRACING_ALLTOALL;
2304 extra->send_size = sendcount;
2305 extra->recv_size = recvcount;
2306 extra->datatype1 = encode_datatype(sendtype);
2307 extra->datatype2 = encode_datatype(recvtype);
2309 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2311 retval = mpi_coll_alltoall_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
2313 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2321 int PMPI_Alltoallv(void *sendbuf, int *sendcounts, int *senddisps,
2322 MPI_Datatype sendtype, void *recvbuf, int *recvcounts,
2323 int *recvdisps, MPI_Datatype recvtype, MPI_Comm comm)
2329 if (comm == MPI_COMM_NULL) {
2330 retval = MPI_ERR_COMM;
2331 } else if (sendtype == MPI_DATATYPE_NULL
2332 || recvtype == MPI_DATATYPE_NULL) {
2333 retval = MPI_ERR_TYPE;
2334 } else if (sendcounts == NULL || senddisps == NULL || recvcounts == NULL
2335 || recvdisps == NULL) {
2336 retval = MPI_ERR_ARG;
2339 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2341 int size = smpi_comm_size(comm);
2342 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2343 extra->type = TRACING_ALLTOALLV;
2344 extra->send_size = 0;
2345 extra->recv_size = 0;
2346 extra->recvcounts= xbt_malloc(size*sizeof(int));
2347 extra->sendcounts= xbt_malloc(size*sizeof(int));
2349 for(i=0; i< size; i++){//copy data to avoid bad free
2350 extra->send_size += sendcounts[i];
2351 extra->recv_size += recvcounts[i];
2353 extra->sendcounts[i] = sendcounts[i];
2354 extra->recvcounts[i] = recvcounts[i];
2356 extra->num_processes = size;
2358 extra->datatype1 = encode_datatype(sendtype);
2359 extra->datatype2 = encode_datatype(recvtype);
2361 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2364 mpi_coll_alltoallv_fun(sendbuf, sendcounts, senddisps, sendtype,
2365 recvbuf, recvcounts, recvdisps, recvtype,
2368 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2377 int PMPI_Get_processor_name(char *name, int *resultlen)
2379 int retval = MPI_SUCCESS;
2381 strncpy(name, SIMIX_host_get_name(SIMIX_host_self()),
2382 strlen(SIMIX_host_get_name(SIMIX_host_self())) < MPI_MAX_PROCESSOR_NAME - 1 ?
2383 strlen(SIMIX_host_get_name(SIMIX_host_self())) +1 :
2384 MPI_MAX_PROCESSOR_NAME - 1 );
2387 MPI_MAX_PROCESSOR_NAME ? MPI_MAX_PROCESSOR_NAME : strlen(name);
2392 int PMPI_Get_count(MPI_Status * status, MPI_Datatype datatype, int *count)
2394 int retval = MPI_SUCCESS;
2397 if (status == NULL || count == NULL) {
2398 retval = MPI_ERR_ARG;
2399 } else if (datatype == MPI_DATATYPE_NULL) {
2400 retval = MPI_ERR_TYPE;
2402 size = smpi_datatype_size(datatype);
2405 } else if (status->count % size != 0) {
2406 retval = MPI_UNDEFINED;
2408 *count = smpi_mpi_get_count(status, datatype);
2414 int PMPI_Type_contiguous(int count, MPI_Datatype old_type, MPI_Datatype* new_type) {
2417 if (old_type == MPI_DATATYPE_NULL) {
2418 retval = MPI_ERR_TYPE;
2419 } else if (count<0){
2420 retval = MPI_ERR_COUNT;
2422 retval = smpi_datatype_contiguous(count, old_type, new_type, 0);
2427 int PMPI_Type_commit(MPI_Datatype* datatype) {
2430 if (datatype == NULL || *datatype == MPI_DATATYPE_NULL) {
2431 retval = MPI_ERR_TYPE;
2433 smpi_datatype_commit(datatype);
2434 retval = MPI_SUCCESS;
2440 int PMPI_Type_vector(int count, int blocklen, int stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2443 if (old_type == MPI_DATATYPE_NULL) {
2444 retval = MPI_ERR_TYPE;
2445 } else if (count<0 || blocklen<0){
2446 retval = MPI_ERR_COUNT;
2448 retval = smpi_datatype_vector(count, blocklen, stride, old_type, new_type);
2453 int PMPI_Type_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2456 if (old_type == MPI_DATATYPE_NULL) {
2457 retval = MPI_ERR_TYPE;
2458 } else if (count<0 || blocklen<0){
2459 retval = MPI_ERR_COUNT;
2461 retval = smpi_datatype_hvector(count, blocklen, stride, old_type, new_type);
2466 int PMPI_Type_create_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2467 return MPI_Type_hvector(count, blocklen, stride, old_type, new_type);
2470 int PMPI_Type_indexed(int count, int* blocklens, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2473 if (old_type == MPI_DATATYPE_NULL) {
2474 retval = MPI_ERR_TYPE;
2475 } else if (count<0){
2476 retval = MPI_ERR_COUNT;
2478 retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
2483 int PMPI_Type_create_indexed(int count, int* blocklens, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2486 if (old_type == MPI_DATATYPE_NULL) {
2487 retval = MPI_ERR_TYPE;
2488 } else if (count<0){
2489 retval = MPI_ERR_COUNT;
2491 retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
2496 int PMPI_Type_create_indexed_block(int count, int blocklength, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2499 if (old_type == MPI_DATATYPE_NULL) {
2500 retval = MPI_ERR_TYPE;
2501 } else if (count<0){
2502 retval = MPI_ERR_COUNT;
2504 int* blocklens=(int*)xbt_malloc(blocklength*count);
2505 for (i=0; i<count;i++)blocklens[i]=blocklength;
2506 retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
2507 xbt_free(blocklens);
2513 int PMPI_Type_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2516 if (old_type == MPI_DATATYPE_NULL) {
2517 retval = MPI_ERR_TYPE;
2518 } else if (count<0){
2519 retval = MPI_ERR_COUNT;
2521 retval = smpi_datatype_hindexed(count, blocklens, indices, old_type, new_type);
2526 int PMPI_Type_create_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2527 return PMPI_Type_hindexed(count, blocklens,indices,old_type,new_type);
2530 int PMPI_Type_create_hindexed_block(int count, int blocklength, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2533 if (old_type == MPI_DATATYPE_NULL) {
2534 retval = MPI_ERR_TYPE;
2535 } else if (count<0){
2536 retval = MPI_ERR_COUNT;
2538 int* blocklens=(int*)xbt_malloc(blocklength*count);
2539 for (i=0; i<count;i++)blocklens[i]=blocklength;
2540 retval = smpi_datatype_hindexed(count, blocklens, indices, old_type, new_type);
2541 xbt_free(blocklens);
2547 int PMPI_Type_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* new_type) {
2551 retval = MPI_ERR_COUNT;
2553 retval = smpi_datatype_struct(count, blocklens, indices, old_types, new_type);
2558 int PMPI_Type_create_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* new_type) {
2559 return PMPI_Type_struct(count, blocklens, indices, old_types, new_type);
2563 int PMPI_Error_class(int errorcode, int* errorclass) {
2564 // assume smpi uses only standard mpi error codes
2565 *errorclass=errorcode;
2570 int PMPI_Initialized(int* flag) {
2571 *flag=smpi_process_initialized();
2575 /* The topo part of MPI_COMM_WORLD should always be NULL. When other topologies
2576 * will be implemented, not only should we check if the topology is NULL, but
2577 * we should check if it is the good topology type (so we have to add a
2578 * MPIR_Topo_Type field, and replace the MPI_Topology field by an union)*/
2580 int PMPI_Cart_create(MPI_Comm comm_old, int ndims, int* dims, int* periodic, int reorder, MPI_Comm* comm_cart) {
2583 if (comm_old == MPI_COMM_NULL){
2584 return MPI_ERR_COMM;
2586 else if (ndims < 0 ||
2587 (ndims > 0 && (dims == NULL ||
2588 periodic == NULL)) ||
2589 comm_cart == NULL) {
2592 retval = smpi_mpi_cart_create(comm_old, ndims, dims, periodic, reorder, comm_cart);
2599 int PMPI_Cart_rank(MPI_Comm comm, int* coords, int* rank) {
2600 if(comm == MPI_COMM_NULL || smpi_comm_topo(comm) == NULL) {
2601 return MPI_ERR_TOPOLOGY;
2603 if (coords == NULL) {
2606 return smpi_mpi_cart_rank(comm, coords, rank);
2609 int PMPI_Cart_shift(MPI_Comm comm, int direction, int displ, int* source, int* dest) {
2610 if(comm == MPI_COMM_NULL || smpi_comm_topo(comm) == NULL) {
2611 return MPI_ERR_TOPOLOGY;
2613 if (source == NULL || dest == NULL || direction < 0 ) {
2616 return smpi_mpi_cart_shift(comm, direction, displ, source, dest);
2619 int PMPI_Cart_coords(MPI_Comm comm, int rank, int maxdims, int* coords) {
2620 if(comm == MPI_COMM_NULL || smpi_comm_topo(comm) == NULL) {
2621 return MPI_ERR_TOPOLOGY;
2623 if (rank < 0 || rank >= smpi_comm_size(comm)) {
2624 return MPI_ERR_RANK;
2629 if(coords == NULL) {
2632 return smpi_mpi_cart_coords(comm, rank, maxdims, coords);
2635 int PMPI_Cart_get(MPI_Comm comm, int maxdims, int* dims, int* periods, int* coords) {
2636 if(comm == NULL || smpi_comm_topo(comm) == NULL) {
2637 return MPI_ERR_TOPOLOGY;
2639 if(maxdims <= 0 || dims == NULL || periods == NULL || coords == NULL) {
2642 return smpi_mpi_cart_get(comm, maxdims, dims, periods, coords);
2645 int PMPI_Cartdim_get(MPI_Comm comm, int* ndims) {
2646 if (comm == MPI_COMM_NULL || smpi_comm_topo(comm) == NULL) {
2647 return MPI_ERR_TOPOLOGY;
2649 if (ndims == NULL) {
2652 return smpi_mpi_cartdim_get(comm, ndims);
2655 int PMPI_Dims_create(int nnodes, int ndims, int* dims) {
2659 if (ndims < 1 || nnodes < 1) {
2660 return MPI_ERR_DIMS;
2663 return smpi_mpi_dims_create(nnodes, ndims, dims);
2666 int PMPI_Cart_sub(MPI_Comm comm, int* remain_dims, MPI_Comm* comm_new) {
2667 if(comm == MPI_COMM_NULL || smpi_comm_topo(comm) == NULL) {
2668 return MPI_ERR_TOPOLOGY;
2670 if (comm_new == NULL) {
2673 return smpi_mpi_cart_sub(comm, remain_dims, comm_new);
2677 /* The following calls are not yet implemented and will fail at runtime. */
2678 /* Once implemented, please move them above this notice. */
2680 #define NOT_YET_IMPLEMENTED { \
2681 XBT_WARN("Not yet implemented : %s. Please contact the Simgrid team if support is needed", __FUNCTION__); \
2682 return MPI_SUCCESS; \
2686 int PMPI_Type_dup(MPI_Datatype datatype, MPI_Datatype *newtype){
2690 int PMPI_Type_set_name(MPI_Datatype datatype, char * name)
2695 int PMPI_Type_get_name(MPI_Datatype datatype, char * name, int* len)
2700 int PMPI_Pack_size(int incount, MPI_Datatype datatype, MPI_Comm comm, int* size) {
2705 int PMPI_Cart_map(MPI_Comm comm_old, int ndims, int* dims, int* periods, int* newrank) {
2710 int PMPI_Graph_create(MPI_Comm comm_old, int nnodes, int* index, int* edges, int reorder, MPI_Comm* comm_graph) {
2714 int PMPI_Graph_get(MPI_Comm comm, int maxindex, int maxedges, int* index, int* edges) {
2718 int PMPI_Graph_map(MPI_Comm comm_old, int nnodes, int* index, int* edges, int* newrank) {
2722 int PMPI_Graph_neighbors(MPI_Comm comm, int rank, int maxneighbors, int* neighbors) {
2726 int PMPI_Graph_neighbors_count(MPI_Comm comm, int rank, int* nneighbors) {
2730 int PMPI_Graphdims_get(MPI_Comm comm, int* nnodes, int* nedges) {
2734 int PMPI_Topo_test(MPI_Comm comm, int* top_type) {
2738 int PMPI_Errhandler_create(MPI_Handler_function* function, MPI_Errhandler* errhandler) {
2742 int PMPI_Errhandler_free(MPI_Errhandler* errhandler) {
2746 int PMPI_Errhandler_get(MPI_Comm comm, MPI_Errhandler* errhandler) {
2750 int PMPI_Error_string(int errorcode, char* string, int* resultlen) {
2754 int PMPI_Errhandler_set(MPI_Comm comm, MPI_Errhandler errhandler) {
2758 int PMPI_Comm_set_errhandler(MPI_Comm comm, MPI_Errhandler errhandler) {
2762 int PMPI_Comm_get_errhandler(MPI_Comm comm, MPI_Errhandler* errhandler) {
2766 int PMPI_Cancel(MPI_Request* request) {
2770 int PMPI_Buffer_attach(void* buffer, int size) {
2774 int PMPI_Buffer_detach(void* buffer, int* size) {
2778 int PMPI_Comm_test_inter(MPI_Comm comm, int* flag) {
2782 int PMPI_Comm_get_attr (MPI_Comm comm, int comm_keyval, void *attribute_val, int *flag)
2787 int PMPI_Comm_set_attr (MPI_Comm comm, int comm_keyval, void *attribute_val)
2792 int PMPI_Comm_delete_attr (MPI_Comm comm, int comm_keyval)
2797 int PMPI_Comm_create_keyval(MPI_Comm_copy_attr_function* copy_fn, MPI_Comm_delete_attr_function* delete_fn, int* keyval, void* extra_state)
2802 int PMPI_Comm_free_keyval(int* keyval) {
2806 int PMPI_Pcontrol(const int level )
2811 int PMPI_Unpack(void* inbuf, int insize, int* position, void* outbuf, int outcount, MPI_Datatype type, MPI_Comm comm) {
2815 int PMPI_Type_get_attr (MPI_Datatype type, int type_keyval, void *attribute_val, int* flag)
2820 int PMPI_Type_set_attr (MPI_Datatype type, int type_keyval, void *attribute_val)
2825 int PMPI_Type_delete_attr (MPI_Datatype type, int comm_keyval)
2830 int PMPI_Type_create_keyval(MPI_Type_copy_attr_function* copy_fn, MPI_Type_delete_attr_function* delete_fn, int* keyval, void* extra_state)
2835 int PMPI_Type_free_keyval(int* keyval) {
2839 int PMPI_Intercomm_create(MPI_Comm local_comm, int local_leader, MPI_Comm peer_comm, int remote_leader, int tag, MPI_Comm* comm_out) {
2843 int PMPI_Intercomm_merge(MPI_Comm comm, int high, MPI_Comm* comm_out) {
2847 int PMPI_Bsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
2851 int PMPI_Bsend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2855 int PMPI_Ibsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2859 int PMPI_Comm_remote_group(MPI_Comm comm, MPI_Group* group) {
2863 int PMPI_Comm_remote_size(MPI_Comm comm, int* size) {
2867 int PMPI_Attr_delete(MPI_Comm comm, int keyval) {
2871 int PMPI_Attr_get(MPI_Comm comm, int keyval, void* attr_value, int* flag) {
2875 int PMPI_Attr_put(MPI_Comm comm, int keyval, void* attr_value) {
2879 int PMPI_Rsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
2883 int PMPI_Rsend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2887 int PMPI_Irsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2891 int PMPI_Keyval_create(MPI_Copy_function* copy_fn, MPI_Delete_function* delete_fn, int* keyval, void* extra_state) {
2895 int PMPI_Keyval_free(int* keyval) {
2899 int PMPI_Test_cancelled(MPI_Status* status, int* flag) {
2903 int PMPI_Pack(void* inbuf, int incount, MPI_Datatype type, void* outbuf, int outcount, int* position, MPI_Comm comm) {
2907 int PMPI_Pack_external_size(char *datarep, int incount, MPI_Datatype datatype, MPI_Aint *size){
2911 int PMPI_Pack_external(char *datarep, void *inbuf, int incount, MPI_Datatype datatype, void *outbuf, MPI_Aint outcount, MPI_Aint *position){
2915 int PMPI_Unpack_external( char *datarep, void *inbuf, MPI_Aint insize, MPI_Aint *position, void *outbuf, int outcount, MPI_Datatype datatype){
2919 int PMPI_Get_elements(MPI_Status* status, MPI_Datatype datatype, int* elements) {
2923 int PMPI_Win_fence( int assert, MPI_Win win){
2927 int PMPI_Win_free( MPI_Win* win){
2931 int PMPI_Win_create( void *base, MPI_Aint size, int disp_unit, MPI_Info info, MPI_Comm comm, MPI_Win *win){
2935 int PMPI_Info_create( MPI_Info *info){
2939 int PMPI_Info_set( MPI_Info info, char *key, char *value){
2943 int PMPI_Info_free( MPI_Info *info){
2947 int PMPI_Get( void *origin_addr, int origin_count, MPI_Datatype origin_datatype, int target_rank,
2948 MPI_Aint target_disp, int target_count, MPI_Datatype target_datatype, MPI_Win win){
2952 int PMPI_Type_get_envelope( MPI_Datatype datatype, int *num_integers,
2953 int *num_addresses, int *num_datatypes, int *combiner){
2957 int PMPI_Type_get_contents(MPI_Datatype datatype, int max_integers, int max_addresses,
2958 int max_datatypes, int* array_of_integers, MPI_Aint* array_of_addresses,
2959 MPI_Datatype* array_of_datatypes){
2963 int PMPI_Type_create_darray(int size, int rank, int ndims, int* array_of_gsizes,
2964 int* array_of_distribs, int* array_of_dargs, int* array_of_psizes,
2965 int order, MPI_Datatype oldtype, MPI_Datatype *newtype) {
2969 int PMPI_Type_create_resized(MPI_Datatype oldtype,MPI_Aint lb, MPI_Aint extent, MPI_Datatype *newtype){
2973 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){
2977 int PMPI_Type_match_size(int typeclass,int size,MPI_Datatype *datatype){
2981 int PMPI_Alltoallw( void *sendbuf, int *sendcnts, int *sdispls, MPI_Datatype *sendtypes,
2982 void *recvbuf, int *recvcnts, int *rdispls, MPI_Datatype *recvtypes,
2987 int PMPI_Comm_set_name(MPI_Comm comm, char* name){
2991 int PMPI_Comm_dup_with_info(MPI_Comm comm, MPI_Info info, MPI_Comm * newcomm){
2995 int PMPI_Comm_split_type(MPI_Comm comm, int split_type, int key, MPI_Info info, MPI_Comm *newcomm){
2999 int PMPI_Comm_set_info (MPI_Comm comm, MPI_Info info){
3003 int PMPI_Comm_get_info (MPI_Comm comm, MPI_Info* info){
3007 int PMPI_Info_get(MPI_Info info,char *key,int valuelen, char *value, int *flag){
3011 int PMPI_Comm_create_errhandler( MPI_Comm_errhandler_fn *function, MPI_Errhandler *errhandler){
3015 int PMPI_Add_error_class( int *errorclass){
3019 int PMPI_Add_error_code( int errorclass, int *errorcode){
3023 int PMPI_Add_error_string( int errorcode, char *string){
3027 int PMPI_Comm_call_errhandler(MPI_Comm comm,int errorcode){
3031 int PMPI_Info_dup(MPI_Info info, MPI_Info *newinfo){
3035 int PMPI_Info_delete(MPI_Info info, char *key){
3039 int PMPI_Info_get_nkeys( MPI_Info info, int *nkeys){
3043 int PMPI_Info_get_nthkey( MPI_Info info, int n, char *key){
3047 int PMPI_Info_get_valuelen( MPI_Info info, char *key, int *valuelen, int *flag){
3051 int PMPI_Request_get_status( MPI_Request request, int *flag, MPI_Status *status){
3055 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){
3059 int PMPI_Grequest_complete( MPI_Request request){
3063 int PMPI_Status_set_cancelled(MPI_Status *status,int flag){
3067 int PMPI_Status_set_elements( MPI_Status *status, MPI_Datatype datatype, int count){
3071 int PMPI_Comm_connect( char *port_name, MPI_Info info, int root, MPI_Comm comm, MPI_Comm *newcomm){
3075 int PMPI_Publish_name( char *service_name, MPI_Info info, char *port_name){
3079 int PMPI_Unpublish_name( char *service_name, MPI_Info info, char *port_name){
3083 int PMPI_Lookup_name( char *service_name, MPI_Info info, char *port_name){
3087 int PMPI_Comm_join( int fd, MPI_Comm *intercomm){
3091 int PMPI_Open_port( MPI_Info info, char *port_name){
3095 int PMPI_Close_port(char *port_name){
3099 int PMPI_Comm_accept( char *port_name, MPI_Info info, int root, MPI_Comm comm, MPI_Comm *newcomm){
3103 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){
3107 int PMPI_Comm_spawn_multiple( int count, char **array_of_commands, char*** array_of_argv,
3108 int* array_of_maxprocs, MPI_Info* array_of_info, int root,
3109 MPI_Comm comm, MPI_Comm *intercomm, int* array_of_errcodes){
3113 int PMPI_Comm_get_parent( MPI_Comm *parent){