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_ARG;
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 smpi_mpi_startall(count, requests);
944 retval = MPI_SUCCESS;
950 int PMPI_Request_free(MPI_Request * request)
955 if (*request == MPI_REQUEST_NULL) {
956 retval = MPI_ERR_ARG;
958 smpi_mpi_request_free(request);
959 retval = MPI_SUCCESS;
965 int PMPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src,
966 int tag, MPI_Comm comm, MPI_Request * request)
972 if (request == NULL) {
973 retval = MPI_ERR_ARG;
974 } else if (comm == MPI_COMM_NULL) {
975 retval = MPI_ERR_COMM;
976 } else if (src == MPI_PROC_NULL) {
977 *request = MPI_REQUEST_NULL;
978 retval = MPI_SUCCESS;
979 } else if (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0)){
980 retval = MPI_ERR_COMM;
981 } else if (count < 0) {
982 retval = MPI_ERR_COUNT;
983 } else if (buf==NULL && count > 0) {
984 retval = MPI_ERR_COUNT;
985 } else if (datatype == MPI_DATATYPE_NULL){
986 retval = MPI_ERR_TYPE;
987 } else if(tag<0 && tag != MPI_ANY_TAG){
988 retval = MPI_ERR_TAG;
992 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
993 int src_traced = smpi_group_index(smpi_comm_group(comm), src);
995 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
996 extra->type = TRACING_IRECV;
997 extra->send_size = count;
998 extra->src = src_traced;
1000 extra->datatype1 = encode_datatype(datatype);
1001 TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__, extra);
1004 *request = smpi_mpi_irecv(buf, count, datatype, src, tag, comm);
1005 retval = MPI_SUCCESS;
1008 TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
1009 (*request)->recv = 1;
1014 if (retval != MPI_SUCCESS && request)
1015 *request = MPI_REQUEST_NULL;
1020 int PMPI_Isend(void *buf, int count, MPI_Datatype datatype, int dst,
1021 int tag, MPI_Comm comm, MPI_Request * request)
1026 if (request == NULL) {
1027 retval = MPI_ERR_ARG;
1028 } else if (comm == MPI_COMM_NULL) {
1029 retval = MPI_ERR_COMM;
1030 } else if (dst == MPI_PROC_NULL) {
1031 *request = MPI_REQUEST_NULL;
1032 retval = MPI_SUCCESS;
1033 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1034 retval = MPI_ERR_RANK;
1035 } else if (count < 0) {
1036 retval = MPI_ERR_COUNT;
1037 } else if (buf==NULL && count > 0) {
1038 retval = MPI_ERR_COUNT;
1039 } else if (datatype == MPI_DATATYPE_NULL){
1040 retval = MPI_ERR_TYPE;
1041 } else if(tag<0 && tag != MPI_ANY_TAG){
1042 retval = MPI_ERR_TAG;
1046 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1047 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1049 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1050 extra->type = TRACING_ISEND;
1051 extra->send_size = count;
1053 extra->dst = dst_traced;
1054 extra->datatype1 = encode_datatype(datatype);
1055 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, extra);
1056 TRACE_smpi_send(rank, rank, dst_traced, count*smpi_datatype_size(datatype));
1059 *request = smpi_mpi_isend(buf, count, datatype, dst, tag, comm);
1060 retval = MPI_SUCCESS;
1063 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1064 (*request)->send = 1;
1069 if (retval != MPI_SUCCESS && request)
1070 *request = MPI_REQUEST_NULL;
1074 int PMPI_Issend(void* buf, int count, MPI_Datatype datatype,
1075 int dst, int tag, MPI_Comm comm, MPI_Request* request)
1080 if (request == NULL) {
1081 retval = MPI_ERR_ARG;
1082 } else if (comm == MPI_COMM_NULL) {
1083 retval = MPI_ERR_COMM;
1084 } else if (dst == MPI_PROC_NULL) {
1085 *request = MPI_REQUEST_NULL;
1086 retval = MPI_SUCCESS;
1087 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1088 retval = MPI_ERR_RANK;
1089 } else if (count < 0) {
1090 retval = MPI_ERR_COUNT;
1091 } else if (buf==NULL && count > 0) {
1092 retval = MPI_ERR_COUNT;
1093 } else if (datatype == MPI_DATATYPE_NULL){
1094 retval = MPI_ERR_TYPE;
1095 } else if(tag<0 && tag != MPI_ANY_TAG){
1096 retval = MPI_ERR_TAG;
1100 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1101 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1102 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1103 extra->type = TRACING_ISSEND;
1104 extra->send_size = count;
1106 extra->dst = dst_traced;
1107 extra->datatype1 = encode_datatype(datatype);
1108 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, extra);
1109 TRACE_smpi_send(rank, rank, dst_traced, count*smpi_datatype_size(datatype));
1112 *request = smpi_mpi_issend(buf, count, datatype, dst, tag, comm);
1113 retval = MPI_SUCCESS;
1116 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1117 (*request)->send = 1;
1122 if (retval != MPI_SUCCESS && request)
1123 *request = MPI_REQUEST_NULL;
1127 int PMPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag,
1128 MPI_Comm comm, MPI_Status * status)
1133 if (comm == MPI_COMM_NULL) {
1134 retval = MPI_ERR_COMM;
1135 } else if (src == MPI_PROC_NULL) {
1136 smpi_empty_status(status);
1137 status->MPI_SOURCE = MPI_PROC_NULL;
1138 retval = MPI_SUCCESS;
1139 } else if (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0)){
1140 retval = MPI_ERR_RANK;
1141 } else if (count < 0) {
1142 retval = MPI_ERR_COUNT;
1143 } else if (buf==NULL && count > 0) {
1144 retval = MPI_ERR_COUNT;
1145 } else if (datatype == MPI_DATATYPE_NULL){
1146 retval = MPI_ERR_TYPE;
1147 } else if(tag<0 && tag != MPI_ANY_TAG){
1148 retval = MPI_ERR_TAG;
1151 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1152 int src_traced = smpi_group_index(smpi_comm_group(comm), src);
1153 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1154 extra->type = TRACING_RECV;
1155 extra->send_size = count;
1156 extra->src = src_traced;
1158 extra->datatype1 = encode_datatype(datatype);
1159 TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__, extra);
1162 smpi_mpi_recv(buf, count, datatype, src, tag, comm, status);
1163 retval = MPI_SUCCESS;
1166 //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1167 if(status!=MPI_STATUS_IGNORE){
1168 src_traced = smpi_group_index(smpi_comm_group(comm), status->MPI_SOURCE);
1169 TRACE_smpi_recv(rank, src_traced, rank);
1171 TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
1179 int PMPI_Send(void *buf, int count, MPI_Datatype datatype, int dst, int tag,
1186 if (comm == MPI_COMM_NULL) {
1187 retval = MPI_ERR_COMM;
1188 } else if (dst == MPI_PROC_NULL) {
1189 retval = MPI_SUCCESS;
1190 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1191 retval = MPI_ERR_RANK;
1192 } else if (count < 0) {
1193 retval = MPI_ERR_COUNT;
1194 } else if (buf==NULL && count > 0) {
1195 retval = MPI_ERR_COUNT;
1196 } else if (datatype == MPI_DATATYPE_NULL){
1197 retval = MPI_ERR_TYPE;
1198 } else if(tag<0 && tag != MPI_ANY_TAG){
1199 retval = MPI_ERR_TAG;
1203 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1204 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1205 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1206 extra->type = TRACING_SEND;
1207 extra->send_size = count;
1209 extra->dst = dst_traced;
1210 extra->datatype1 = encode_datatype(datatype);
1211 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, extra);
1212 TRACE_smpi_send(rank, rank, dst_traced,count*smpi_datatype_size(datatype));
1215 smpi_mpi_send(buf, count, datatype, dst, tag, comm);
1216 retval = MPI_SUCCESS;
1219 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1229 int PMPI_Ssend(void* buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm) {
1234 if (comm == MPI_COMM_NULL) {
1235 retval = MPI_ERR_COMM;
1236 } else if (dst == MPI_PROC_NULL) {
1237 retval = MPI_SUCCESS;
1238 } else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0){
1239 retval = MPI_ERR_RANK;
1240 } else if (count < 0) {
1241 retval = MPI_ERR_COUNT;
1242 } else if (buf==NULL && count > 0) {
1243 retval = MPI_ERR_COUNT;
1244 } else if (datatype == MPI_DATATYPE_NULL){
1245 retval = MPI_ERR_TYPE;
1246 } else if(tag<0 && tag != MPI_ANY_TAG){
1247 retval = MPI_ERR_TAG;
1251 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1252 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1253 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1254 extra->type = TRACING_SSEND;
1255 extra->send_size = count;
1257 extra->dst = dst_traced;
1258 extra->datatype1 = encode_datatype(datatype);
1259 TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, extra); TRACE_smpi_send(rank, rank, dst_traced,count*smpi_datatype_size(datatype));
1262 smpi_mpi_ssend(buf, count, datatype, dst, tag, comm);
1263 retval = MPI_SUCCESS;
1266 TRACE_smpi_ptp_out(rank, rank, dst_traced, __FUNCTION__);
1274 int PMPI_Sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1275 int dst, int sendtag, void *recvbuf, int recvcount,
1276 MPI_Datatype recvtype, int src, int recvtag,
1277 MPI_Comm comm, MPI_Status * status)
1283 if (comm == MPI_COMM_NULL) {
1284 retval = MPI_ERR_COMM;
1285 } else if (sendtype == MPI_DATATYPE_NULL
1286 || recvtype == MPI_DATATYPE_NULL) {
1287 retval = MPI_ERR_TYPE;
1288 } else if (src == MPI_PROC_NULL || dst == MPI_PROC_NULL) {
1289 smpi_empty_status(status);
1290 status->MPI_SOURCE = MPI_PROC_NULL;
1291 retval = MPI_SUCCESS;
1292 }else if (dst >= smpi_group_size(smpi_comm_group(comm)) || dst <0 ||
1293 (src!=MPI_ANY_SOURCE && (src >= smpi_group_size(smpi_comm_group(comm)) || src <0))){
1294 retval = MPI_ERR_RANK;
1295 } else if (sendcount < 0 || recvcount<0) {
1296 retval = MPI_ERR_COUNT;
1297 } else if ((sendbuf==NULL && sendcount > 0)||(recvbuf==NULL && recvcount>0)) {
1298 retval = MPI_ERR_COUNT;
1299 } else if((sendtag<0 && sendtag != MPI_ANY_TAG)||(recvtag<0 && recvtag != MPI_ANY_TAG)){
1300 retval = MPI_ERR_TAG;
1304 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1305 int dst_traced = smpi_group_index(smpi_comm_group(comm), dst);
1306 int src_traced = smpi_group_index(smpi_comm_group(comm), src);
1307 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1308 extra->type = TRACING_SENDRECV;
1309 extra->send_size = sendcount;
1310 extra->recv_size = recvcount;
1311 extra->src = src_traced;
1312 extra->dst = dst_traced;
1313 extra->datatype1 = encode_datatype(sendtype);
1314 extra->datatype2 = encode_datatype(recvtype);
1316 TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__, extra);
1317 TRACE_smpi_send(rank, rank, dst_traced,sendcount*smpi_datatype_size(sendtype));
1321 smpi_mpi_sendrecv(sendbuf, sendcount, sendtype, dst, sendtag, recvbuf,
1322 recvcount, recvtype, src, recvtag, comm, status);
1323 retval = MPI_SUCCESS;
1326 TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1327 TRACE_smpi_recv(rank, src_traced, rank);
1336 int PMPI_Sendrecv_replace(void *buf, int count, MPI_Datatype datatype,
1337 int dst, int sendtag, int src, int recvtag,
1338 MPI_Comm comm, MPI_Status * status)
1340 //TODO: suboptimal implementation
1343 if (datatype == MPI_DATATYPE_NULL) {
1344 retval = MPI_ERR_TYPE;
1345 } else if (count < 0) {
1346 retval = MPI_ERR_COUNT;
1348 int size = smpi_datatype_get_extent(datatype) * count;
1349 recvbuf = xbt_new0(char, size);
1351 MPI_Sendrecv(buf, count, datatype, dst, sendtag, recvbuf, count,
1352 datatype, src, recvtag, comm, status);
1353 if(retval==MPI_SUCCESS){
1354 smpi_datatype_copy(recvbuf, count, datatype, buf, count, datatype);
1362 int PMPI_Test(MPI_Request * request, int *flag, MPI_Status * status)
1366 if (request == NULL || flag == NULL) {
1367 retval = MPI_ERR_ARG;
1368 } else if (*request == MPI_REQUEST_NULL) {
1370 smpi_empty_status(status);
1371 retval = MPI_ERR_REQUEST;
1374 int rank = request && (*request)->comm != MPI_COMM_NULL
1375 ? smpi_process_index()
1378 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1379 extra->type = TRACING_TEST;
1380 TRACE_smpi_testing_in(rank, extra);
1382 *flag = smpi_mpi_test(request, status);
1384 TRACE_smpi_testing_out(rank);
1386 retval = MPI_SUCCESS;
1392 int PMPI_Testany(int count, MPI_Request requests[], int *index, int *flag,
1393 MPI_Status * status)
1398 if (index == NULL || flag == NULL) {
1399 retval = MPI_ERR_ARG;
1401 *flag = smpi_mpi_testany(count, requests, index, status);
1402 retval = MPI_SUCCESS;
1408 int PMPI_Testall(int count, MPI_Request* requests, int* flag, MPI_Status* statuses)
1414 retval = MPI_ERR_ARG;
1416 *flag = smpi_mpi_testall(count, requests, statuses);
1417 retval = MPI_SUCCESS;
1423 int PMPI_Probe(int source, int tag, MPI_Comm comm, MPI_Status* status) {
1427 if (status == NULL) {
1428 retval = MPI_ERR_ARG;
1429 } else if (comm == MPI_COMM_NULL) {
1430 retval = MPI_ERR_COMM;
1431 } else if (source == MPI_PROC_NULL) {
1432 smpi_empty_status(status);
1433 status->MPI_SOURCE = MPI_PROC_NULL;
1434 retval = MPI_SUCCESS;
1436 smpi_mpi_probe(source, tag, comm, status);
1437 retval = MPI_SUCCESS;
1444 int PMPI_Iprobe(int source, int tag, MPI_Comm comm, int* flag, MPI_Status* status) {
1449 retval = MPI_ERR_ARG;
1450 } else if (status == NULL) {
1451 retval = MPI_ERR_ARG;
1452 } else if (comm == MPI_COMM_NULL) {
1453 retval = MPI_ERR_COMM;
1454 } else if (source == MPI_PROC_NULL) {
1456 smpi_empty_status(status);
1457 status->MPI_SOURCE = MPI_PROC_NULL;
1458 retval = MPI_SUCCESS;
1460 smpi_mpi_iprobe(source, tag, comm, flag, status);
1461 retval = MPI_SUCCESS;
1467 int PMPI_Wait(MPI_Request * request, MPI_Status * status)
1473 smpi_empty_status(status);
1475 if (request == NULL) {
1476 retval = MPI_ERR_ARG;
1477 } else if (*request == MPI_REQUEST_NULL) {
1478 retval = MPI_ERR_REQUEST;
1482 int rank = request && (*request)->comm != MPI_COMM_NULL
1483 ? smpi_process_index()
1486 int src_traced = (*request)->src;
1487 int dst_traced = (*request)->dst;
1488 MPI_Comm comm = (*request)->comm;
1489 int is_wait_for_receive = (*request)->recv;
1490 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1491 extra->type = TRACING_WAIT;
1492 TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__, extra);
1495 smpi_mpi_wait(request, status);
1496 retval = MPI_SUCCESS;
1499 //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1500 TRACE_smpi_ptp_out(rank, src_traced, dst_traced, __FUNCTION__);
1501 if (is_wait_for_receive) {
1502 if(src_traced==MPI_ANY_SOURCE)
1503 src_traced = (status!=MPI_STATUS_IGNORE) ?
1504 smpi_group_rank(smpi_comm_group(comm), status->MPI_SOURCE) :
1506 TRACE_smpi_recv(rank, src_traced, dst_traced);
1516 int PMPI_Waitany(int count, MPI_Request requests[], int *index, MPI_Status * status)
1523 //save requests information for tracing
1525 int *srcs = xbt_new0(int, count);
1526 int *dsts = xbt_new0(int, count);
1527 int *recvs = xbt_new0(int, count);
1528 MPI_Comm *comms = xbt_new0(MPI_Comm, count);
1530 for (i = 0; i < count; i++) {
1531 MPI_Request req = requests[i]; //already received requests are no longer valid
1535 recvs[i] = req->recv;
1536 comms[i] = req->comm;
1539 int rank_traced = smpi_process_index();
1540 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1541 extra->type = TRACING_WAITANY;
1542 extra->send_size=count;
1543 TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__,extra);
1546 *index = smpi_mpi_waitany(count, requests, status);
1548 if(*index!=MPI_UNDEFINED){
1549 int src_traced = srcs[*index];
1550 //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1551 int dst_traced = dsts[*index];
1552 int is_wait_for_receive = recvs[*index];
1553 if (is_wait_for_receive) {
1554 if(srcs[*index]==MPI_ANY_SOURCE)
1555 src_traced = (status!=MPI_STATUSES_IGNORE) ?
1556 smpi_group_rank(smpi_comm_group(comms[*index]), status->MPI_SOURCE) :
1558 TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1560 TRACE_smpi_ptp_out(rank_traced, src_traced, dst_traced, __FUNCTION__);
1572 int PMPI_Waitall(int count, MPI_Request requests[], MPI_Status status[])
1577 //save information from requests
1579 int *srcs = xbt_new0(int, count);
1580 int *dsts = xbt_new0(int, count);
1581 int *recvs = xbt_new0(int, count);
1582 int *valid = xbt_new0(int, count);
1583 MPI_Comm *comms = xbt_new0(MPI_Comm, count);
1585 //int valid_count = 0;
1586 for (i = 0; i < count; i++) {
1587 MPI_Request req = requests[i];
1588 if(req!=MPI_REQUEST_NULL){
1591 recvs[i] = req->recv;
1592 comms[i] = req->comm;
1598 int rank_traced = smpi_process_index();
1599 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1600 extra->type = TRACING_WAITALL;
1601 extra->send_size=count;
1602 TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__,extra);
1604 int retval = smpi_mpi_waitall(count, requests, status);
1606 for (i = 0; i < count; i++) {
1608 //int src_traced = srcs[*index];
1609 //the src may not have been known at the beginning of the recv (MPI_ANY_SOURCE)
1610 int src_traced = srcs[i];
1611 int dst_traced = dsts[i];
1612 int is_wait_for_receive = recvs[i];
1613 if (is_wait_for_receive) {
1614 if(src_traced==MPI_ANY_SOURCE)
1615 src_traced = (status!=MPI_STATUSES_IGNORE) ?
1616 smpi_group_rank(smpi_comm_group(comms[i]), status[i].MPI_SOURCE) :
1618 TRACE_smpi_recv(rank_traced, src_traced, dst_traced);
1622 TRACE_smpi_ptp_out(rank_traced, -1, -1, __FUNCTION__);
1634 int PMPI_Waitsome(int incount, MPI_Request requests[], int *outcount,
1635 int *indices, MPI_Status status[])
1640 if (outcount == NULL) {
1641 retval = MPI_ERR_ARG;
1643 *outcount = smpi_mpi_waitsome(incount, requests, indices, status);
1644 retval = MPI_SUCCESS;
1650 int PMPI_Testsome(int incount, MPI_Request requests[], int* outcount,
1651 int* indices, MPI_Status status[])
1656 if (outcount == NULL) {
1657 retval = MPI_ERR_ARG;
1659 *outcount = smpi_mpi_testsome(incount, requests, indices, status);
1660 retval = MPI_SUCCESS;
1667 int PMPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root, MPI_Comm comm)
1673 if (comm == MPI_COMM_NULL) {
1674 retval = MPI_ERR_COMM;
1677 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1678 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1680 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1681 extra->type = TRACING_BCAST;
1682 extra->send_size = count;
1683 extra->root = root_traced;
1684 extra->datatype1 = encode_datatype(datatype);
1685 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__, extra);
1688 mpi_coll_bcast_fun(buf, count, datatype, root, comm);
1689 retval = MPI_SUCCESS;
1691 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1699 int PMPI_Barrier(MPI_Comm comm)
1705 if (comm == MPI_COMM_NULL) {
1706 retval = MPI_ERR_COMM;
1709 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1710 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1711 extra->type = TRACING_BARRIER;
1712 TRACE_smpi_collective_in(rank, -1, __FUNCTION__, extra);
1714 mpi_coll_barrier_fun(comm);
1715 retval = MPI_SUCCESS;
1717 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1725 int PMPI_Gather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1726 void *recvbuf, int recvcount, MPI_Datatype recvtype,
1727 int root, MPI_Comm comm)
1733 if (comm == MPI_COMM_NULL) {
1734 retval = MPI_ERR_COMM;
1735 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1736 ((smpi_comm_rank(comm) == root) && (recvtype == MPI_DATATYPE_NULL))){
1737 retval = MPI_ERR_TYPE;
1738 } else if ((( sendbuf != MPI_IN_PLACE) && (sendcount <0)) ||
1739 ((smpi_comm_rank(comm) == root) && (recvcount <0))){
1740 retval = MPI_ERR_COUNT;
1743 char* sendtmpbuf = (char*) sendbuf;
1744 int sendtmpcount = sendcount;
1745 MPI_Datatype sendtmptype = sendtype;
1746 if( (smpi_comm_rank(comm) == root) && (sendbuf == MPI_IN_PLACE )) {
1748 sendtmptype=recvtype;
1751 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1752 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1753 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1754 extra->type = TRACING_GATHER;
1755 extra->send_size = sendtmpcount;
1756 extra->recv_size = recvcount;
1757 extra->root = root_traced;
1758 extra->datatype1 = encode_datatype(sendtmptype);
1759 extra->datatype2 = encode_datatype(recvtype);
1761 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__, extra);
1763 mpi_coll_gather_fun(sendtmpbuf, sendtmpcount, sendtmptype, recvbuf, recvcount,
1764 recvtype, root, comm);
1767 retval = MPI_SUCCESS;
1769 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1777 int PMPI_Gatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1778 void *recvbuf, int *recvcounts, int *displs,
1779 MPI_Datatype recvtype, int root, MPI_Comm comm)
1785 if (comm == MPI_COMM_NULL) {
1786 retval = MPI_ERR_COMM;
1787 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1788 ((smpi_comm_rank(comm) == root) && (recvtype == MPI_DATATYPE_NULL))){
1789 retval = MPI_ERR_TYPE;
1790 } else if (( sendbuf != MPI_IN_PLACE) && (sendcount <0)){
1791 retval = MPI_ERR_COUNT;
1792 } else if (recvcounts == NULL || displs == NULL) {
1793 retval = MPI_ERR_ARG;
1795 char* sendtmpbuf = (char*) sendbuf;
1796 int sendtmpcount = sendcount;
1797 MPI_Datatype sendtmptype = sendtype;
1798 if( (smpi_comm_rank(comm) == root) && (sendbuf == MPI_IN_PLACE )) {
1800 sendtmptype=recvtype;
1804 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1805 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1807 int size = smpi_comm_size(comm);
1808 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1809 extra->type = TRACING_GATHERV;
1810 extra->send_size = sendtmpcount;
1811 extra->recvcounts= xbt_malloc(size*sizeof(int));
1812 for(i=0; i< size; i++)//copy data to avoid bad free
1813 extra->recvcounts[i] = recvcounts[i];
1814 extra->num_processes = size;
1815 extra->root = root_traced;
1816 extra->datatype1 = encode_datatype(sendtmptype);
1817 extra->datatype2 = encode_datatype(recvtype);
1819 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__,extra);
1821 smpi_mpi_gatherv(sendtmpbuf, sendtmpcount, sendtmptype, recvbuf, recvcounts,
1822 displs, recvtype, root, comm);
1823 retval = MPI_SUCCESS;
1825 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1833 int PMPI_Allgather(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1834 void *recvbuf, int recvcount, MPI_Datatype recvtype,
1841 if (comm == MPI_COMM_NULL) {
1842 retval = MPI_ERR_COMM;
1843 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1844 (recvtype == MPI_DATATYPE_NULL)){
1845 retval = MPI_ERR_TYPE;
1846 } else if ((( sendbuf != MPI_IN_PLACE) && (sendcount <0)) ||
1848 retval = MPI_ERR_COUNT;
1850 if(sendbuf == MPI_IN_PLACE) {
1851 sendbuf=((char*)recvbuf)+smpi_datatype_get_extent(recvtype)*recvcount*smpi_comm_rank(comm);
1852 sendcount=recvcount;
1856 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1857 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1858 extra->type = TRACING_ALLGATHER;
1859 extra->send_size = sendcount;
1860 extra->recv_size = recvcount;
1861 extra->datatype1 = encode_datatype(sendtype);
1862 extra->datatype2 = encode_datatype(recvtype);
1864 TRACE_smpi_collective_in(rank, -1, __FUNCTION__, extra);
1866 mpi_coll_allgather_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1868 retval = MPI_SUCCESS;
1871 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1878 int PMPI_Allgatherv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1879 void *recvbuf, int *recvcounts, int *displs,
1880 MPI_Datatype recvtype, MPI_Comm comm)
1886 if (comm == MPI_COMM_NULL) {
1887 retval = MPI_ERR_COMM;
1888 } else if ((( sendbuf != MPI_IN_PLACE) && (sendtype == MPI_DATATYPE_NULL)) ||
1889 (recvtype == MPI_DATATYPE_NULL)){
1890 retval = MPI_ERR_TYPE;
1891 } else if (( sendbuf != MPI_IN_PLACE) && (sendcount <0)){
1892 retval = MPI_ERR_COUNT;
1893 } else if (recvcounts == NULL || displs == NULL) {
1894 retval = MPI_ERR_ARG;
1897 if(sendbuf == MPI_IN_PLACE) {
1898 sendbuf=((char*)recvbuf)+smpi_datatype_get_extent(recvtype)*displs[smpi_comm_rank(comm)];
1899 sendcount=recvcounts[smpi_comm_rank(comm)];
1903 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1905 int size = smpi_comm_size(comm);
1906 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1907 extra->type = TRACING_ALLGATHERV;
1908 extra->send_size = sendcount;
1909 extra->recvcounts= xbt_malloc(size*sizeof(int));
1910 for(i=0; i< size; i++)//copy data to avoid bad free
1911 extra->recvcounts[i] = recvcounts[i];
1912 extra->num_processes = size;
1913 extra->datatype1 = encode_datatype(sendtype);
1914 extra->datatype2 = encode_datatype(recvtype);
1916 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
1918 mpi_coll_allgatherv_fun(sendbuf, sendcount, sendtype, recvbuf, recvcounts,
1919 displs, recvtype, comm);
1920 retval = MPI_SUCCESS;
1922 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
1930 int PMPI_Scatter(void *sendbuf, int sendcount, MPI_Datatype sendtype,
1931 void *recvbuf, int recvcount, MPI_Datatype recvtype,
1932 int root, MPI_Comm comm)
1938 if (comm == MPI_COMM_NULL) {
1939 retval = MPI_ERR_COMM;
1940 } else if (((smpi_comm_rank(comm)==root) && (sendtype == MPI_DATATYPE_NULL))
1941 || ((recvbuf !=MPI_IN_PLACE) && (recvtype == MPI_DATATYPE_NULL))) {
1942 retval = MPI_ERR_TYPE;
1945 if (recvbuf == MPI_IN_PLACE) {
1947 recvcount=sendcount;
1950 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1951 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1952 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
1953 extra->type = TRACING_SCATTER;
1954 extra->send_size = sendcount;
1955 extra->recv_size= recvcount;
1956 extra->root = root_traced;
1957 extra->datatype1 = encode_datatype(sendtype);
1958 extra->datatype2 = encode_datatype(recvtype);
1960 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__,extra);
1962 mpi_coll_scatter_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount,
1963 recvtype, root, comm);
1964 retval = MPI_SUCCESS;
1966 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
1974 int PMPI_Scatterv(void *sendbuf, int *sendcounts, int *displs,
1975 MPI_Datatype sendtype, void *recvbuf, int recvcount,
1976 MPI_Datatype recvtype, int root, MPI_Comm comm)
1982 if (comm == MPI_COMM_NULL) {
1983 retval = MPI_ERR_COMM;
1984 } else if (sendcounts == NULL || displs == NULL) {
1985 retval = MPI_ERR_ARG;
1986 } else if (((smpi_comm_rank(comm)==root) && (sendtype == MPI_DATATYPE_NULL))
1987 || ((recvbuf !=MPI_IN_PLACE) && (recvtype == MPI_DATATYPE_NULL))) {
1988 retval = MPI_ERR_TYPE;
1990 if (recvbuf == MPI_IN_PLACE) {
1992 recvcount=sendcounts[smpi_comm_rank(comm)];
1995 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
1996 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
1998 int size = smpi_comm_size(comm);
1999 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2000 extra->type = TRACING_SCATTERV;
2001 extra->recv_size = recvcount;
2002 extra->sendcounts= xbt_malloc(size*sizeof(int));
2003 for(i=0; i< size; i++)//copy data to avoid bad free
2004 extra->sendcounts[i] = sendcounts[i];
2005 extra->num_processes = size;
2006 extra->root = root_traced;
2007 extra->datatype1 = encode_datatype(sendtype);
2008 extra->datatype2 = encode_datatype(recvtype);
2010 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__,extra);
2013 smpi_mpi_scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf,
2014 recvcount, recvtype, root, comm);
2015 retval = MPI_SUCCESS;
2017 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
2025 int PMPI_Reduce(void *sendbuf, void *recvbuf, int count,
2026 MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm)
2032 if (comm == MPI_COMM_NULL) {
2033 retval = MPI_ERR_COMM;
2034 } else if (datatype == MPI_DATATYPE_NULL || op == MPI_OP_NULL) {
2035 retval = MPI_ERR_ARG;
2038 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2039 int root_traced = smpi_group_index(smpi_comm_group(comm), root);
2040 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2041 extra->type = TRACING_REDUCE;
2042 extra->send_size = count;
2043 extra->datatype1 = encode_datatype(datatype);
2044 extra->root = root_traced;
2046 TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__,extra);
2048 mpi_coll_reduce_fun(sendbuf, recvbuf, count, datatype, op, root, comm);
2050 retval = MPI_SUCCESS;
2052 TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__);
2060 int PMPI_Reduce_local(void *inbuf, void *inoutbuf, int count,
2061 MPI_Datatype datatype, MPI_Op op){
2065 if (datatype == MPI_DATATYPE_NULL || op == MPI_OP_NULL) {
2066 retval = MPI_ERR_ARG;
2068 smpi_op_apply(op, inbuf, inoutbuf, &count, &datatype);
2075 int PMPI_Allreduce(void *sendbuf, void *recvbuf, int count,
2076 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2082 if (comm == MPI_COMM_NULL) {
2083 retval = MPI_ERR_COMM;
2084 } else if (datatype == MPI_DATATYPE_NULL) {
2085 retval = MPI_ERR_TYPE;
2086 } else if (op == MPI_OP_NULL) {
2087 retval = MPI_ERR_OP;
2090 char* sendtmpbuf = (char*) sendbuf;
2091 if( sendbuf == MPI_IN_PLACE ) {
2092 sendtmpbuf = (char *)xbt_malloc(count*smpi_datatype_get_extent(datatype));
2093 smpi_datatype_copy(recvbuf, count, datatype,sendtmpbuf, count, datatype);
2096 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2097 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2098 extra->type = TRACING_ALLREDUCE;
2099 extra->send_size = count;
2100 extra->datatype1 = encode_datatype(datatype);
2102 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2104 mpi_coll_allreduce_fun(sendtmpbuf, recvbuf, count, datatype, op, comm);
2106 if( sendbuf == MPI_IN_PLACE ) {
2107 xbt_free(sendtmpbuf);
2110 retval = MPI_SUCCESS;
2112 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2120 int PMPI_Scan(void *sendbuf, void *recvbuf, int count,
2121 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2127 if (comm == MPI_COMM_NULL) {
2128 retval = MPI_ERR_COMM;
2129 } else if (datatype == MPI_DATATYPE_NULL) {
2130 retval = MPI_ERR_TYPE;
2131 } else if (op == MPI_OP_NULL) {
2132 retval = MPI_ERR_OP;
2135 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2136 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2137 extra->type = TRACING_SCAN;
2138 extra->send_size = count;
2139 extra->datatype1 = encode_datatype(datatype);
2141 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2143 smpi_mpi_scan(sendbuf, recvbuf, count, datatype, op, comm);
2144 retval = MPI_SUCCESS;
2146 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2154 int PMPI_Exscan(void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype,
2155 MPI_Op op, MPI_Comm comm){
2160 if (comm == MPI_COMM_NULL) {
2161 retval = MPI_ERR_COMM;
2162 } else if (datatype == MPI_DATATYPE_NULL) {
2163 retval = MPI_ERR_TYPE;
2164 } else if (op == MPI_OP_NULL) {
2165 retval = MPI_ERR_OP;
2168 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2169 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2170 extra->type = TRACING_EXSCAN;
2171 extra->send_size = count;
2172 extra->datatype1 = encode_datatype(datatype);
2174 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2176 smpi_mpi_exscan(sendbuf, recvbuf, count, datatype, op, comm);
2177 retval = MPI_SUCCESS;
2179 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2187 int PMPI_Reduce_scatter(void *sendbuf, void *recvbuf, int *recvcounts,
2188 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2193 if (comm == MPI_COMM_NULL) {
2194 retval = MPI_ERR_COMM;
2195 } else if (datatype == MPI_DATATYPE_NULL) {
2196 retval = MPI_ERR_TYPE;
2197 } else if (op == MPI_OP_NULL) {
2198 retval = MPI_ERR_OP;
2199 } else if (recvcounts == NULL) {
2200 retval = MPI_ERR_ARG;
2203 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2205 int size = smpi_comm_size(comm);
2206 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2207 extra->type = TRACING_REDUCE_SCATTER;
2208 extra->send_size = 0;
2209 extra->recvcounts= xbt_malloc(size*sizeof(int));
2210 for(i=0; i< size; i++)//copy data to avoid bad free
2211 extra->recvcounts[i] = recvcounts[i];
2212 extra->num_processes = size;
2213 extra->datatype1 = encode_datatype(datatype);
2215 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2217 void* sendtmpbuf=sendbuf;
2218 if(sendbuf==MPI_IN_PLACE){
2222 mpi_coll_reduce_scatter_fun(sendtmpbuf, recvbuf, recvcounts,
2223 datatype, op, comm);
2224 retval = MPI_SUCCESS;
2226 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2234 int PMPI_Reduce_scatter_block(void *sendbuf, void *recvbuf, int recvcount,
2235 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
2240 if (comm == MPI_COMM_NULL) {
2241 retval = MPI_ERR_COMM;
2242 } else if (datatype == MPI_DATATYPE_NULL) {
2243 retval = MPI_ERR_TYPE;
2244 } else if (op == MPI_OP_NULL) {
2245 retval = MPI_ERR_OP;
2246 } else if (recvcount < 0) {
2247 retval = MPI_ERR_ARG;
2249 int count=smpi_comm_size(comm);
2252 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2253 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2254 extra->type = TRACING_REDUCE_SCATTER;
2255 extra->send_size = 0;
2256 extra->recvcounts= xbt_malloc(count*sizeof(int));
2257 for(i=0; i< count; i++)//copy data to avoid bad free
2258 extra->recvcounts[i] = recvcount;
2259 extra->num_processes = count;
2260 extra->datatype1 = encode_datatype(datatype);
2262 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2264 int* recvcounts=(int*)xbt_malloc(count);
2265 for (i=0; i<count;i++)recvcounts[i]=recvcount;
2266 mpi_coll_reduce_scatter_fun(sendbuf, recvbuf, recvcounts,
2267 datatype, op, comm);
2268 xbt_free(recvcounts);
2269 retval = MPI_SUCCESS;
2271 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2279 int PMPI_Alltoall(void *sendbuf, int sendcount, MPI_Datatype sendtype,
2280 void *recvbuf, int recvcount, MPI_Datatype recvtype,
2287 if (comm == MPI_COMM_NULL) {
2288 retval = MPI_ERR_COMM;
2289 } else if (sendtype == MPI_DATATYPE_NULL
2290 || recvtype == MPI_DATATYPE_NULL) {
2291 retval = MPI_ERR_TYPE;
2294 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2295 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2296 extra->type = TRACING_ALLTOALL;
2297 extra->send_size = sendcount;
2298 extra->recv_size = recvcount;
2299 extra->datatype1 = encode_datatype(sendtype);
2300 extra->datatype2 = encode_datatype(recvtype);
2302 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2304 retval = mpi_coll_alltoall_fun(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm);
2306 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2314 int PMPI_Alltoallv(void *sendbuf, int *sendcounts, int *senddisps,
2315 MPI_Datatype sendtype, void *recvbuf, int *recvcounts,
2316 int *recvdisps, MPI_Datatype recvtype, MPI_Comm comm)
2322 if (comm == MPI_COMM_NULL) {
2323 retval = MPI_ERR_COMM;
2324 } else if (sendtype == MPI_DATATYPE_NULL
2325 || recvtype == MPI_DATATYPE_NULL) {
2326 retval = MPI_ERR_TYPE;
2327 } else if (sendcounts == NULL || senddisps == NULL || recvcounts == NULL
2328 || recvdisps == NULL) {
2329 retval = MPI_ERR_ARG;
2332 int rank = comm != MPI_COMM_NULL ? smpi_process_index() : -1;
2334 int size = smpi_comm_size(comm);
2335 instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
2336 extra->type = TRACING_ALLTOALLV;
2337 extra->send_size = 0;
2338 extra->recv_size = 0;
2339 extra->recvcounts= xbt_malloc(size*sizeof(int));
2340 extra->sendcounts= xbt_malloc(size*sizeof(int));
2342 for(i=0; i< size; i++){//copy data to avoid bad free
2343 extra->send_size += sendcounts[i];
2344 extra->recv_size += recvcounts[i];
2346 extra->sendcounts[i] = sendcounts[i];
2347 extra->recvcounts[i] = recvcounts[i];
2349 extra->num_processes = size;
2351 extra->datatype1 = encode_datatype(sendtype);
2352 extra->datatype2 = encode_datatype(recvtype);
2354 TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
2357 mpi_coll_alltoallv_fun(sendbuf, sendcounts, senddisps, sendtype,
2358 recvbuf, recvcounts, recvdisps, recvtype,
2361 TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
2370 int PMPI_Get_processor_name(char *name, int *resultlen)
2372 int retval = MPI_SUCCESS;
2374 strncpy(name, SIMIX_host_get_name(SIMIX_host_self()),
2375 strlen(SIMIX_host_get_name(SIMIX_host_self())) < MPI_MAX_PROCESSOR_NAME - 1 ?
2376 strlen(SIMIX_host_get_name(SIMIX_host_self())) +1 :
2377 MPI_MAX_PROCESSOR_NAME - 1 );
2380 MPI_MAX_PROCESSOR_NAME ? MPI_MAX_PROCESSOR_NAME : strlen(name);
2385 int PMPI_Get_count(MPI_Status * status, MPI_Datatype datatype, int *count)
2387 int retval = MPI_SUCCESS;
2390 if (status == NULL || count == NULL) {
2391 retval = MPI_ERR_ARG;
2392 } else if (datatype == MPI_DATATYPE_NULL) {
2393 retval = MPI_ERR_TYPE;
2395 size = smpi_datatype_size(datatype);
2398 } else if (status->count % size != 0) {
2399 retval = MPI_UNDEFINED;
2401 *count = smpi_mpi_get_count(status, datatype);
2407 int PMPI_Type_contiguous(int count, MPI_Datatype old_type, MPI_Datatype* new_type) {
2410 if (old_type == MPI_DATATYPE_NULL) {
2411 retval = MPI_ERR_TYPE;
2412 } else if (count<0){
2413 retval = MPI_ERR_COUNT;
2415 retval = smpi_datatype_contiguous(count, old_type, new_type, 0);
2420 int PMPI_Type_commit(MPI_Datatype* datatype) {
2423 if (datatype == NULL || *datatype == MPI_DATATYPE_NULL) {
2424 retval = MPI_ERR_TYPE;
2426 smpi_datatype_commit(datatype);
2427 retval = MPI_SUCCESS;
2433 int PMPI_Type_vector(int count, int blocklen, int stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2436 if (old_type == MPI_DATATYPE_NULL) {
2437 retval = MPI_ERR_TYPE;
2438 } else if (count<0 || blocklen<0){
2439 retval = MPI_ERR_COUNT;
2441 retval = smpi_datatype_vector(count, blocklen, stride, old_type, new_type);
2446 int PMPI_Type_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2449 if (old_type == MPI_DATATYPE_NULL) {
2450 retval = MPI_ERR_TYPE;
2451 } else if (count<0 || blocklen<0){
2452 retval = MPI_ERR_COUNT;
2454 retval = smpi_datatype_hvector(count, blocklen, stride, old_type, new_type);
2459 int PMPI_Type_create_hvector(int count, int blocklen, MPI_Aint stride, MPI_Datatype old_type, MPI_Datatype* new_type) {
2460 return MPI_Type_hvector(count, blocklen, stride, old_type, new_type);
2463 int PMPI_Type_indexed(int count, int* blocklens, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2466 if (old_type == MPI_DATATYPE_NULL) {
2467 retval = MPI_ERR_TYPE;
2468 } else if (count<0){
2469 retval = MPI_ERR_COUNT;
2471 retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
2476 int PMPI_Type_create_indexed(int count, int* blocklens, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2479 if (old_type == MPI_DATATYPE_NULL) {
2480 retval = MPI_ERR_TYPE;
2481 } else if (count<0){
2482 retval = MPI_ERR_COUNT;
2484 retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
2489 int PMPI_Type_create_indexed_block(int count, int blocklength, int* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2492 if (old_type == MPI_DATATYPE_NULL) {
2493 retval = MPI_ERR_TYPE;
2494 } else if (count<0){
2495 retval = MPI_ERR_COUNT;
2497 int* blocklens=(int*)xbt_malloc(blocklength*count);
2498 for (i=0; i<count;i++)blocklens[i]=blocklength;
2499 retval = smpi_datatype_indexed(count, blocklens, indices, old_type, new_type);
2500 xbt_free(blocklens);
2506 int PMPI_Type_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2509 if (old_type == MPI_DATATYPE_NULL) {
2510 retval = MPI_ERR_TYPE;
2511 } else if (count<0){
2512 retval = MPI_ERR_COUNT;
2514 retval = smpi_datatype_hindexed(count, blocklens, indices, old_type, new_type);
2519 int PMPI_Type_create_hindexed(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2520 return PMPI_Type_hindexed(count, blocklens,indices,old_type,new_type);
2523 int PMPI_Type_create_hindexed_block(int count, int blocklength, MPI_Aint* indices, MPI_Datatype old_type, MPI_Datatype* new_type) {
2526 if (old_type == MPI_DATATYPE_NULL) {
2527 retval = MPI_ERR_TYPE;
2528 } else if (count<0){
2529 retval = MPI_ERR_COUNT;
2531 int* blocklens=(int*)xbt_malloc(blocklength*count);
2532 for (i=0; i<count;i++)blocklens[i]=blocklength;
2533 retval = smpi_datatype_hindexed(count, blocklens, indices, old_type, new_type);
2534 xbt_free(blocklens);
2540 int PMPI_Type_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* new_type) {
2544 retval = MPI_ERR_COUNT;
2546 retval = smpi_datatype_struct(count, blocklens, indices, old_types, new_type);
2551 int PMPI_Type_create_struct(int count, int* blocklens, MPI_Aint* indices, MPI_Datatype* old_types, MPI_Datatype* new_type) {
2552 return PMPI_Type_struct(count, blocklens, indices, old_types, new_type);
2556 int PMPI_Error_class(int errorcode, int* errorclass) {
2557 // assume smpi uses only standard mpi error codes
2558 *errorclass=errorcode;
2563 int PMPI_Initialized(int* flag) {
2564 *flag=smpi_process_initialized();
2568 /* The topo part of MPI_COMM_WORLD should always be NULL. When other topologies
2569 * will be implemented, not only should we check if the topology is NULL, but
2570 * we should check if it is the good topology type (so we have to add a
2571 * MPIR_Topo_Type field, and replace the MPI_Topology field by an union)*/
2573 int PMPI_Cart_create(MPI_Comm comm_old, int ndims, int* dims, int* periodic, int reorder, MPI_Comm* comm_cart) {
2576 if (comm_old == MPI_COMM_NULL){
2577 return MPI_ERR_COMM;
2579 else if (ndims < 0 ||
2580 (ndims > 0 && (dims == NULL ||
2581 periodic == NULL)) ||
2582 comm_cart == NULL) {
2585 retval = smpi_mpi_cart_create(comm_old, ndims, dims, periodic, reorder, comm_cart);
2592 int PMPI_Cart_rank(MPI_Comm comm, int* coords, int* rank) {
2593 if(comm == MPI_COMM_NULL || smpi_comm_topo(comm) == NULL) {
2594 return MPI_ERR_TOPOLOGY;
2596 if (coords == NULL) {
2599 return smpi_mpi_cart_rank(comm, coords, rank);
2602 int PMPI_Cart_shift(MPI_Comm comm, int direction, int displ, int* source, int* dest) {
2603 if(comm == MPI_COMM_NULL || smpi_comm_topo(comm) == NULL) {
2604 return MPI_ERR_TOPOLOGY;
2606 if (source == NULL || dest == NULL || direction < 0 ) {
2609 return smpi_mpi_cart_shift(comm, direction, displ, source, dest);
2612 int PMPI_Cart_coords(MPI_Comm comm, int rank, int maxdims, int* coords) {
2613 if(comm == MPI_COMM_NULL || smpi_comm_topo(comm) == NULL) {
2614 return MPI_ERR_TOPOLOGY;
2616 if (rank < 0 || rank >= smpi_comm_size(comm)) {
2617 return MPI_ERR_RANK;
2622 if(coords == NULL) {
2625 return smpi_mpi_cart_coords(comm, rank, maxdims, coords);
2628 int PMPI_Cart_get(MPI_Comm comm, int maxdims, int* dims, int* periods, int* coords) {
2629 if(comm == NULL || smpi_comm_topo(comm) == NULL) {
2630 return MPI_ERR_TOPOLOGY;
2632 if(maxdims <= 0 || dims == NULL || periods == NULL || coords == NULL) {
2635 return smpi_mpi_cart_get(comm, maxdims, dims, periods, coords);
2638 int PMPI_Cartdim_get(MPI_Comm comm, int* ndims) {
2639 if (comm == MPI_COMM_NULL || smpi_comm_topo(comm) == NULL) {
2640 return MPI_ERR_TOPOLOGY;
2642 if (ndims == NULL) {
2645 return smpi_mpi_cartdim_get(comm, ndims);
2648 int PMPI_Dims_create(int nnodes, int ndims, int* dims) {
2652 if (ndims < 1 || nnodes < 1) {
2653 return MPI_ERR_DIMS;
2656 return smpi_mpi_dims_create(nnodes, ndims, dims);
2659 int PMPI_Cart_sub(MPI_Comm comm, int* remain_dims, MPI_Comm* comm_new) {
2660 if(comm == MPI_COMM_NULL || smpi_comm_topo(comm) == NULL) {
2661 return MPI_ERR_TOPOLOGY;
2663 if (comm_new == NULL) {
2666 return smpi_mpi_cart_sub(comm, remain_dims, comm_new);
2670 /* The following calls are not yet implemented and will fail at runtime. */
2671 /* Once implemented, please move them above this notice. */
2673 #define NOT_YET_IMPLEMENTED { \
2674 XBT_WARN("Not yet implemented : %s. Please contact the Simgrid team if support is needed", __FUNCTION__); \
2675 return MPI_SUCCESS; \
2679 int PMPI_Type_dup(MPI_Datatype datatype, MPI_Datatype *newtype){
2683 int PMPI_Type_set_name(MPI_Datatype datatype, char * name)
2688 int PMPI_Type_get_name(MPI_Datatype datatype, char * name, int* len)
2693 int PMPI_Pack_size(int incount, MPI_Datatype datatype, MPI_Comm comm, int* size) {
2698 int PMPI_Cart_map(MPI_Comm comm_old, int ndims, int* dims, int* periods, int* newrank) {
2703 int PMPI_Graph_create(MPI_Comm comm_old, int nnodes, int* index, int* edges, int reorder, MPI_Comm* comm_graph) {
2707 int PMPI_Graph_get(MPI_Comm comm, int maxindex, int maxedges, int* index, int* edges) {
2711 int PMPI_Graph_map(MPI_Comm comm_old, int nnodes, int* index, int* edges, int* newrank) {
2715 int PMPI_Graph_neighbors(MPI_Comm comm, int rank, int maxneighbors, int* neighbors) {
2719 int PMPI_Graph_neighbors_count(MPI_Comm comm, int rank, int* nneighbors) {
2723 int PMPI_Graphdims_get(MPI_Comm comm, int* nnodes, int* nedges) {
2727 int PMPI_Topo_test(MPI_Comm comm, int* top_type) {
2731 int PMPI_Errhandler_create(MPI_Handler_function* function, MPI_Errhandler* errhandler) {
2735 int PMPI_Errhandler_free(MPI_Errhandler* errhandler) {
2739 int PMPI_Errhandler_get(MPI_Comm comm, MPI_Errhandler* errhandler) {
2743 int PMPI_Error_string(int errorcode, char* string, int* resultlen) {
2747 int PMPI_Errhandler_set(MPI_Comm comm, MPI_Errhandler errhandler) {
2751 int PMPI_Comm_set_errhandler(MPI_Comm comm, MPI_Errhandler errhandler) {
2755 int PMPI_Comm_get_errhandler(MPI_Comm comm, MPI_Errhandler* errhandler) {
2759 int PMPI_Cancel(MPI_Request* request) {
2763 int PMPI_Buffer_attach(void* buffer, int size) {
2767 int PMPI_Buffer_detach(void* buffer, int* size) {
2771 int PMPI_Comm_test_inter(MPI_Comm comm, int* flag) {
2775 int PMPI_Comm_get_attr (MPI_Comm comm, int comm_keyval, void *attribute_val, int *flag)
2780 int PMPI_Comm_set_attr (MPI_Comm comm, int comm_keyval, void *attribute_val)
2785 int PMPI_Comm_delete_attr (MPI_Comm comm, int comm_keyval)
2790 int PMPI_Comm_create_keyval(MPI_Comm_copy_attr_function* copy_fn, MPI_Comm_delete_attr_function* delete_fn, int* keyval, void* extra_state)
2795 int PMPI_Comm_free_keyval(int* keyval) {
2799 int PMPI_Pcontrol(const int level )
2804 int PMPI_Unpack(void* inbuf, int insize, int* position, void* outbuf, int outcount, MPI_Datatype type, MPI_Comm comm) {
2808 int PMPI_Type_get_attr (MPI_Datatype type, int type_keyval, void *attribute_val, int* flag)
2813 int PMPI_Type_set_attr (MPI_Datatype type, int type_keyval, void *attribute_val)
2818 int PMPI_Type_delete_attr (MPI_Datatype type, int comm_keyval)
2823 int PMPI_Type_create_keyval(MPI_Type_copy_attr_function* copy_fn, MPI_Type_delete_attr_function* delete_fn, int* keyval, void* extra_state)
2828 int PMPI_Type_free_keyval(int* keyval) {
2832 int PMPI_Intercomm_create(MPI_Comm local_comm, int local_leader, MPI_Comm peer_comm, int remote_leader, int tag, MPI_Comm* comm_out) {
2836 int PMPI_Intercomm_merge(MPI_Comm comm, int high, MPI_Comm* comm_out) {
2840 int PMPI_Bsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
2844 int PMPI_Bsend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2848 int PMPI_Ibsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2852 int PMPI_Comm_remote_group(MPI_Comm comm, MPI_Group* group) {
2856 int PMPI_Comm_remote_size(MPI_Comm comm, int* size) {
2860 int PMPI_Attr_delete(MPI_Comm comm, int keyval) {
2864 int PMPI_Attr_get(MPI_Comm comm, int keyval, void* attr_value, int* flag) {
2868 int PMPI_Attr_put(MPI_Comm comm, int keyval, void* attr_value) {
2872 int PMPI_Rsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm) {
2876 int PMPI_Rsend_init(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2880 int PMPI_Irsend(void* buf, int count, MPI_Datatype datatype, int dest, int tag, MPI_Comm comm, MPI_Request* request) {
2884 int PMPI_Keyval_create(MPI_Copy_function* copy_fn, MPI_Delete_function* delete_fn, int* keyval, void* extra_state) {
2888 int PMPI_Keyval_free(int* keyval) {
2892 int PMPI_Test_cancelled(MPI_Status* status, int* flag) {
2896 int PMPI_Pack(void* inbuf, int incount, MPI_Datatype type, void* outbuf, int outcount, int* position, MPI_Comm comm) {
2900 int PMPI_Pack_external_size(char *datarep, int incount, MPI_Datatype datatype, MPI_Aint *size){
2904 int PMPI_Pack_external(char *datarep, void *inbuf, int incount, MPI_Datatype datatype, void *outbuf, MPI_Aint outcount, MPI_Aint *position){
2908 int PMPI_Unpack_external( char *datarep, void *inbuf, MPI_Aint insize, MPI_Aint *position, void *outbuf, int outcount, MPI_Datatype datatype){
2912 int PMPI_Get_elements(MPI_Status* status, MPI_Datatype datatype, int* elements) {
2916 int PMPI_Win_fence( int assert, MPI_Win win){
2920 int PMPI_Win_free( MPI_Win* win){
2924 int PMPI_Win_create( void *base, MPI_Aint size, int disp_unit, MPI_Info info, MPI_Comm comm, MPI_Win *win){
2928 int PMPI_Info_create( MPI_Info *info){
2932 int PMPI_Info_set( MPI_Info info, char *key, char *value){
2936 int PMPI_Info_free( MPI_Info *info){
2940 int PMPI_Get( void *origin_addr, int origin_count, MPI_Datatype origin_datatype, int target_rank,
2941 MPI_Aint target_disp, int target_count, MPI_Datatype target_datatype, MPI_Win win){
2945 int PMPI_Type_get_envelope( MPI_Datatype datatype, int *num_integers,
2946 int *num_addresses, int *num_datatypes, int *combiner){
2950 int PMPI_Type_get_contents(MPI_Datatype datatype, int max_integers, int max_addresses,
2951 int max_datatypes, int* array_of_integers, MPI_Aint* array_of_addresses,
2952 MPI_Datatype* array_of_datatypes){
2956 int PMPI_Type_create_darray(int size, int rank, int ndims, int* array_of_gsizes,
2957 int* array_of_distribs, int* array_of_dargs, int* array_of_psizes,
2958 int order, MPI_Datatype oldtype, MPI_Datatype *newtype) {
2962 int PMPI_Type_create_resized(MPI_Datatype oldtype,MPI_Aint lb, MPI_Aint extent, MPI_Datatype *newtype){
2966 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){
2970 int PMPI_Type_match_size(int typeclass,int size,MPI_Datatype *datatype){
2974 int PMPI_Alltoallw( void *sendbuf, int *sendcnts, int *sdispls, MPI_Datatype *sendtypes,
2975 void *recvbuf, int *recvcnts, int *rdispls, MPI_Datatype *recvtypes,
2980 int PMPI_Comm_set_name(MPI_Comm comm, char* name){
2984 int PMPI_Comm_dup_with_info(MPI_Comm comm, MPI_Info info, MPI_Comm * newcomm){
2988 int PMPI_Comm_split_type(MPI_Comm comm, int split_type, int key, MPI_Info info, MPI_Comm *newcomm){
2992 int PMPI_Comm_set_info (MPI_Comm comm, MPI_Info info){
2996 int PMPI_Comm_get_info (MPI_Comm comm, MPI_Info* info){
3000 int PMPI_Info_get(MPI_Info info,char *key,int valuelen, char *value, int *flag){
3004 int PMPI_Comm_create_errhandler( MPI_Comm_errhandler_fn *function, MPI_Errhandler *errhandler){
3008 int PMPI_Add_error_class( int *errorclass){
3012 int PMPI_Add_error_code( int errorclass, int *errorcode){
3016 int PMPI_Add_error_string( int errorcode, char *string){
3020 int PMPI_Comm_call_errhandler(MPI_Comm comm,int errorcode){
3024 int PMPI_Info_dup(MPI_Info info, MPI_Info *newinfo){
3028 int PMPI_Info_delete(MPI_Info info, char *key){
3032 int PMPI_Info_get_nkeys( MPI_Info info, int *nkeys){
3036 int PMPI_Info_get_nthkey( MPI_Info info, int n, char *key){
3040 int PMPI_Info_get_valuelen( MPI_Info info, char *key, int *valuelen, int *flag){
3044 int PMPI_Request_get_status( MPI_Request request, int *flag, MPI_Status *status){
3048 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){
3052 int PMPI_Grequest_complete( MPI_Request request){
3056 int PMPI_Status_set_cancelled(MPI_Status *status,int flag){
3060 int PMPI_Status_set_elements( MPI_Status *status, MPI_Datatype datatype, int count){
3064 int PMPI_Comm_connect( char *port_name, MPI_Info info, int root, MPI_Comm comm, MPI_Comm *newcomm){
3068 int PMPI_Publish_name( char *service_name, MPI_Info info, char *port_name){
3072 int PMPI_Unpublish_name( char *service_name, MPI_Info info, char *port_name){
3076 int PMPI_Lookup_name( char *service_name, MPI_Info info, char *port_name){
3080 int PMPI_Comm_join( int fd, MPI_Comm *intercomm){
3084 int PMPI_Open_port( MPI_Info info, char *port_name){
3088 int PMPI_Close_port(char *port_name){
3092 int PMPI_Comm_accept( char *port_name, MPI_Info info, int root, MPI_Comm comm, MPI_Comm *newcomm){
3096 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){
3100 int PMPI_Comm_spawn_multiple( int count, char **array_of_commands, char*** array_of_argv,
3101 int* array_of_maxprocs, MPI_Info* array_of_info, int root,
3102 MPI_Comm comm, MPI_Comm *intercomm, int* array_of_errcodes){
3106 int PMPI_Comm_get_parent( MPI_Comm *parent){