X-Git-Url: http://info.iut-bm.univ-fcomte.fr/pub/gitweb/simgrid.git/blobdiff_plain/04444a3077f46c717719a05b0c10b39d16a8c53a..b2bb0172d45aeb7a2e873ce982e6174279772b91:/src/smpi/smpi_replay.c diff --git a/src/smpi/smpi_replay.c b/src/smpi/smpi_replay.c index 3d2f1be123..9151825c05 100644 --- a/src/smpi/smpi_replay.c +++ b/src/smpi/smpi_replay.c @@ -77,6 +77,33 @@ static MPI_Datatype decode_datatype(const char *const action) return MPI_CURRENT_TYPE; } + +const char* encode_datatype(MPI_Datatype datatype) +{ + + //default type for output is set to MPI_BYTE + // MPI_DEFAULT_TYPE is not set for output, use directly MPI_BYTE + if (datatype==MPI_BYTE){ + return ""; + } + if(datatype==MPI_DOUBLE) + return "0"; + if(datatype==MPI_INT) + return "1"; + if(datatype==MPI_CHAR) + return "2"; + if(datatype==MPI_SHORT) + return "3"; + if(datatype==MPI_LONG) + return "4"; + if(datatype==MPI_FLOAT) + return "5"; + + // default - not implemented. + // do not warn here as we pass in this function even for other trace formats + return "-1"; +} + static void action_init(const char *const *action) { int i; @@ -140,7 +167,18 @@ static void action_comm_dup(const char *const *action) static void action_compute(const char *const *action) { double clock = smpi_process_simulated_elapsed(); - smpi_execute_flops(parse_double(action[2])); + double flops= parse_double(action[2]); +#ifdef HAVE_TRACING + int rank = smpi_comm_rank(MPI_COMM_WORLD); + instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1); + extra->type=TRACING_COMPUTING; + extra->comp_size=flops; + TRACE_smpi_computing_in(rank, extra); +#endif + smpi_execute_flops(flops); +#ifdef HAVE_TRACING + TRACE_smpi_computing_out(rank); +#endif log_timed_action (action, clock); } @@ -161,7 +199,13 @@ static void action_send(const char *const *action) int rank = smpi_comm_rank(MPI_COMM_WORLD); int dst_traced = smpi_group_rank(smpi_comm_group(MPI_COMM_WORLD), to); - TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, size*smpi_datatype_size(MPI_CURRENT_TYPE)); + instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1); + extra->type = TRACING_SEND; + extra->send_size = size; + extra->src = rank; + extra->dst = dst_traced; + extra->datatype1 = encode_datatype(MPI_CURRENT_TYPE); + TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, extra); TRACE_smpi_send(rank, rank, dst_traced, size*smpi_datatype_size(MPI_CURRENT_TYPE)); #endif @@ -188,7 +232,13 @@ static void action_Isend(const char *const *action) #ifdef HAVE_TRACING int rank = smpi_comm_rank(MPI_COMM_WORLD); int dst_traced = smpi_group_rank(smpi_comm_group(MPI_COMM_WORLD), to); - TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, size*smpi_datatype_size(MPI_CURRENT_TYPE)); + instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1); + extra->type = TRACING_ISEND; + extra->send_size = size; + extra->src = rank; + extra->dst = dst_traced; + extra->datatype1 = encode_datatype(MPI_CURRENT_TYPE); + TRACE_smpi_ptp_in(rank, rank, dst_traced, __FUNCTION__, extra); TRACE_smpi_send(rank, rank, dst_traced, size*smpi_datatype_size(MPI_CURRENT_TYPE)); #endif @@ -217,7 +267,13 @@ static void action_recv(const char *const *action) { int rank = smpi_comm_rank(MPI_COMM_WORLD); int src_traced = smpi_group_rank(smpi_comm_group(MPI_COMM_WORLD), from); - TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__, size*smpi_datatype_size(MPI_CURRENT_TYPE)); + instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1); + extra->type = TRACING_RECV; + extra->send_size = size; + extra->src = src_traced; + extra->dst = rank; + extra->datatype1 = encode_datatype(MPI_CURRENT_TYPE); + TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__, extra); #endif smpi_mpi_recv(NULL, size, MPI_CURRENT_TYPE, from, 0, MPI_COMM_WORLD, &status); @@ -246,7 +302,13 @@ static void action_Irecv(const char *const *action) #ifdef HAVE_TRACING int rank = smpi_comm_rank(MPI_COMM_WORLD); int src_traced = smpi_group_rank(smpi_comm_group(MPI_COMM_WORLD), from); - TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__, size*smpi_datatype_size(MPI_CURRENT_TYPE)); + instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1); + extra->type = TRACING_IRECV; + extra->send_size = size; + extra->src = src_traced; + extra->dst = rank; + extra->datatype1 = encode_datatype(MPI_CURRENT_TYPE); + TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__, extra); #endif request = smpi_mpi_irecv(NULL, size, MPI_CURRENT_TYPE, from, 0, MPI_COMM_WORLD); @@ -281,7 +343,9 @@ static void action_wait(const char *const *action){ int src_traced = smpi_group_rank(group, request->src); int dst_traced = smpi_group_rank(group, request->dst); int is_wait_for_receive = request->recv; - TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__, -1); + instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1); + extra->type = TRACING_WAIT; + TRACE_smpi_ptp_in(rank, src_traced, dst_traced, __FUNCTION__, extra); #endif smpi_mpi_wait(&request, &status); #ifdef HAVE_TRACING @@ -338,8 +402,11 @@ static void action_waitall(const char *const *action){ } } int rank_traced = smpi_process_index(); - TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__, count_requests); - #endif + instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1); + extra->type = TRACING_WAITALL; + extra->send_size=count_requests; + TRACE_smpi_ptp_in(rank_traced, -1, -1, __FUNCTION__,extra); + #endif smpi_mpi_waitall(count_requests, requests, status); @@ -369,7 +436,9 @@ static void action_barrier(const char *const *action){ double clock = smpi_process_simulated_elapsed(); #ifdef HAVE_TRACING int rank = smpi_comm_rank(MPI_COMM_WORLD); - TRACE_smpi_collective_in(rank, -1, __FUNCTION__, smpi_comm_size(MPI_COMM_WORLD)); + instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1); + extra->type = TRACING_BARRIER; + TRACE_smpi_collective_in(rank, -1, __FUNCTION__, extra); #endif smpi_mpi_barrier(MPI_COMM_WORLD); #ifdef HAVE_TRACING @@ -400,11 +469,18 @@ static void action_bcast(const char *const *action) #ifdef HAVE_TRACING int rank = smpi_comm_rank(MPI_COMM_WORLD); - int root_traced = smpi_group_rank(smpi_comm_group(MPI_COMM_WORLD), 0); - TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__,size*smpi_datatype_size(MPI_CURRENT_TYPE)); + int root_traced = smpi_group_index(smpi_comm_group(MPI_COMM_WORLD), root); + + instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1); + extra->type = TRACING_BCAST; + extra->send_size = size; + extra->root = root_traced; + extra->datatype1 = encode_datatype(MPI_CURRENT_TYPE); + TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__, extra); + #endif - smpi_mpi_bcast(NULL, size, MPI_CURRENT_TYPE, root, MPI_COMM_WORLD); + mpi_coll_bcast_fun(NULL, size, MPI_CURRENT_TYPE, root, MPI_COMM_WORLD); #ifdef HAVE_TRACING TRACE_smpi_collective_out(rank, root_traced, __FUNCTION__); #endif @@ -429,8 +505,15 @@ static void action_reduce(const char *const *action) #ifdef HAVE_TRACING int rank = smpi_comm_rank(MPI_COMM_WORLD); - int root_traced = smpi_group_rank(smpi_comm_group(MPI_COMM_WORLD), 0); - TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__,comm_size*smpi_datatype_size(MPI_CURRENT_TYPE)); + int root_traced = smpi_group_rank(smpi_comm_group(MPI_COMM_WORLD), root); + instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1); + extra->type = TRACING_REDUCE; + extra->send_size = comm_size; + extra->comp_size = comp_size; + extra->datatype1 = encode_datatype(MPI_CURRENT_TYPE); + extra->root = root_traced; + + TRACE_smpi_collective_in(rank, root_traced, __FUNCTION__,extra); #endif mpi_coll_reduce_fun(NULL, NULL, comm_size, MPI_CURRENT_TYPE, MPI_OP_NULL, root, MPI_COMM_WORLD); smpi_execute_flops(comp_size); @@ -451,7 +534,13 @@ static void action_allReduce(const char *const *action) { double clock = smpi_process_simulated_elapsed(); #ifdef HAVE_TRACING int rank = smpi_comm_rank(MPI_COMM_WORLD); - TRACE_smpi_collective_in(rank, -1, __FUNCTION__,comp_size*smpi_datatype_size(MPI_CURRENT_TYPE)); + instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1); + extra->type = TRACING_ALLREDUCE; + extra->send_size = comm_size; + extra->comp_size = comp_size; + extra->datatype1 = encode_datatype(MPI_CURRENT_TYPE); + + TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra); #endif mpi_coll_reduce_fun(NULL, NULL, comm_size, MPI_CURRENT_TYPE, MPI_OP_NULL, 0, MPI_COMM_WORLD); smpi_execute_flops(comp_size); @@ -483,7 +572,14 @@ static void action_allToAll(const char *const *action) { #ifdef HAVE_TRACING int rank = smpi_process_index(); - TRACE_smpi_collective_in(rank, -1, __FUNCTION__,send_size*smpi_datatype_size(MPI_CURRENT_TYPE)); + instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1); + extra->type = TRACING_ALLTOALL; + extra->send_size = send_size; + extra->recv_size = recv_size; + extra->datatype1 = encode_datatype(MPI_CURRENT_TYPE); + extra->datatype2 = encode_datatype(MPI_CURRENT_TYPE2); + + TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra); #endif mpi_coll_alltoall_fun(send, send_size, MPI_CURRENT_TYPE, recv, recv_size, MPI_CURRENT_TYPE2, MPI_COMM_WORLD); @@ -524,7 +620,7 @@ static void action_gather(const char *const *action) { MPI_CURRENT_TYPE2=MPI_DEFAULT_TYPE; } void *send = calloc(send_size, smpi_datatype_size(MPI_CURRENT_TYPE)); - void *recv = calloc(recv_size, smpi_datatype_size(MPI_CURRENT_TYPE2)); + void *recv = NULL; int root=atoi(action[4]); int rank = smpi_process_index(); @@ -533,7 +629,15 @@ static void action_gather(const char *const *action) { recv = calloc(recv_size*comm_size, smpi_datatype_size(MPI_CURRENT_TYPE2)); #ifdef HAVE_TRACING - TRACE_smpi_collective_in(rank, -1, __FUNCTION__,send_size*smpi_datatype_size(MPI_CURRENT_TYPE)); + instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1); + extra->type = TRACING_GATHER; + extra->send_size = send_size; + extra->recv_size = recv_size; + extra->root = root; + extra->datatype1 = encode_datatype(MPI_CURRENT_TYPE); + extra->datatype2 = encode_datatype(MPI_CURRENT_TYPE2); + + TRACE_smpi_collective_in(rank, root, __FUNCTION__, extra); #endif smpi_mpi_gather(send, send_size, MPI_CURRENT_TYPE, recv, recv_size, MPI_CURRENT_TYPE2, @@ -549,6 +653,79 @@ smpi_mpi_gather(send, send_size, MPI_CURRENT_TYPE, } + +static void action_gatherv(const char *const *action) { + /* + The structure of the gatherv action for the rank 0 (total 4 processes) + is the following: + 0 gather 68 68 10 10 10 0 0 0 + + where: + 1) 68 is the sendcount + 2) 68 10 10 10 is the recvcounts + 3) 0 is the root node + 4) 0 is the send datatype id, see decode_datatype() + 5) 0 is the recv datatype id, see decode_datatype() + */ + double clock = smpi_process_simulated_elapsed(); + int comm_size = smpi_comm_size(MPI_COMM_WORLD); + int send_size = parse_double(action[2]); + int *disps = xbt_new0(int, comm_size); + int *recvcounts = xbt_new0(int, comm_size); + int i=0,recv_sum=0; + + MPI_Datatype MPI_CURRENT_TYPE2; + if(action[4+comm_size]) { + MPI_CURRENT_TYPE=decode_datatype(action[4+comm_size]); + MPI_CURRENT_TYPE2=decode_datatype(action[5+comm_size]); + } else { + MPI_CURRENT_TYPE=MPI_DEFAULT_TYPE; + MPI_CURRENT_TYPE2=MPI_DEFAULT_TYPE; + } + void *send = calloc(send_size, smpi_datatype_size(MPI_CURRENT_TYPE)); + void *recv = NULL; + for(i=0;itype = TRACING_GATHERV; + extra->send_size = send_size; + extra->recvcounts= xbt_malloc(comm_size*sizeof(int)); + for(i=0; i< comm_size; i++)//copy data to avoid bad free + extra->recvcounts[i] = recvcounts[i]; + extra->root = root; + extra->num_processes = comm_size; + extra->datatype1 = encode_datatype(MPI_CURRENT_TYPE); + extra->datatype2 = encode_datatype(MPI_CURRENT_TYPE2); + + TRACE_smpi_collective_in(rank, root, __FUNCTION__, extra); +#endif +smpi_mpi_gatherv(send, send_size, MPI_CURRENT_TYPE, + recv, recvcounts, disps, MPI_CURRENT_TYPE2, + root, MPI_COMM_WORLD); + +#ifdef HAVE_TRACING + TRACE_smpi_collective_out(rank, -1, __FUNCTION__); +#endif + + log_timed_action (action, clock); + xbt_free(recvcounts); + xbt_free(send); + xbt_free(recv); + xbt_free(disps); + +} + static void action_reducescatter(const char *const *action) { /* @@ -586,7 +763,18 @@ static void action_reducescatter(const char *const *action) { } #ifdef HAVE_TRACING - TRACE_smpi_collective_in(rank, -1, __FUNCTION__, recv_sum*smpi_datatype_size(MPI_CURRENT_TYPE)); + instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1); + extra->type = TRACING_REDUCE_SCATTER; + extra->send_size = 0; + extra->recvcounts= xbt_malloc(comm_size*sizeof(int)); + for(i=0; i< comm_size; i++)//copy data to avoid bad free + extra->recvcounts[i] = recvcounts[i]; + extra->datatype1 = encode_datatype(MPI_CURRENT_TYPE); + extra->comp_size = comp_size; + extra->num_processes = comm_size; + + + TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra); #endif mpi_coll_reduce_fun(NULL, NULL, recv_sum, MPI_CURRENT_TYPE, MPI_OP_NULL, root, MPI_COMM_WORLD); @@ -598,7 +786,8 @@ static void action_reducescatter(const char *const *action) { #ifdef HAVE_TRACING TRACE_smpi_collective_out(rank, -1, __FUNCTION__); #endif - + xbt_free(recvcounts); + xbt_free(disps); log_timed_action (action, clock); } @@ -644,8 +833,18 @@ static void action_allgatherv(const char *const *action) { void *recvbuf = calloc(recv_sum, smpi_datatype_size(MPI_CURRENT_TYPE2)); #ifdef HAVE_TRACING - int rank = MPI_COMM_WORLD != MPI_COMM_NULL ? smpi_process_index() : -1; - TRACE_smpi_collective_in(rank, -1, __FUNCTION__,sendcount*smpi_datatype_size(MPI_CURRENT_TYPE)); + int rank = smpi_process_index(); + instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1); + extra->type = TRACING_ALLGATHERV; + extra->send_size = sendcount; + extra->recvcounts= xbt_malloc(comm_size*sizeof(int)); + for(i=0; i< comm_size; i++)//copy data to avoid bad free + extra->recvcounts[i] = recvcounts[i]; + extra->datatype1 = encode_datatype(MPI_CURRENT_TYPE); + extra->datatype2 = encode_datatype(MPI_CURRENT_TYPE2); + extra->num_processes = comm_size; + + TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra); #endif mpi_coll_allgatherv_fun(sendbuf, sendcount, MPI_CURRENT_TYPE, recvbuf, recvcounts, disps, MPI_CURRENT_TYPE2, MPI_COMM_WORLD); @@ -709,10 +908,23 @@ static void action_allToAllv(const char *const *action) { #ifdef HAVE_TRACING - int rank = MPI_COMM_WORLD != MPI_COMM_NULL ? smpi_process_index() : -1; - int count=0; - for(i=0;itype = TRACING_ALLTOALLV; + extra->recvcounts= xbt_malloc(comm_size*sizeof(int)); + extra->sendcounts= xbt_malloc(comm_size*sizeof(int)); + extra->num_processes = comm_size; + + for(i=0; i< comm_size; i++){//copy data to avoid bad free + extra->send_size += sendcounts[i]; + extra->sendcounts[i] = sendcounts[i]; + extra->recv_size += recvcounts[i]; + extra->recvcounts[i] = recvcounts[i]; + } + extra->datatype1 = encode_datatype(MPI_CURRENT_TYPE); + extra->datatype2 = encode_datatype(MPI_CURRENT_TYPE2); + + TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra); #endif mpi_coll_alltoallv_fun(sendbuf, sendcounts, senddisps, MPI_CURRENT_TYPE, recvbuf, recvcounts, recvdisps, MPI_CURRENT_TYPE, @@ -731,7 +943,18 @@ static void action_allToAllv(const char *const *action) { } void smpi_replay_init(int *argc, char***argv){ - PMPI_Init(argc, argv); + smpi_process_init(argc, argv); + smpi_process_mark_as_initialized(); +#ifdef HAVE_TRACING + int rank = smpi_process_index(); + TRACE_smpi_init(rank); + TRACE_smpi_computing_init(rank); + instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1); + extra->type = TRACING_INIT; + TRACE_smpi_collective_in(rank, -1, __FUNCTION__, extra); + TRACE_smpi_collective_out(rank, -1, __FUNCTION__); +#endif + if (!smpi_process_index()){ _xbt_replay_action_init(); xbt_replay_action_register("init", action_init); @@ -752,6 +975,7 @@ void smpi_replay_init(int *argc, char***argv){ xbt_replay_action_register("allToAll", action_allToAll); xbt_replay_action_register("allToAllV", action_allToAllv); xbt_replay_action_register("gather", action_gather); + xbt_replay_action_register("gatherV", action_gatherv); xbt_replay_action_register("allGatherV", action_allgatherv); xbt_replay_action_register("reduceScatter", action_reducescatter); xbt_replay_action_register("compute", action_compute); @@ -777,5 +1001,17 @@ int smpi_replay_finalize(){ reqq = NULL; } smpi_mpi_barrier(MPI_COMM_WORLD); - return PMPI_Finalize(); +#ifdef HAVE_TRACING + int rank = smpi_process_index(); + instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1); + extra->type = TRACING_FINALIZE; + TRACE_smpi_collective_in(rank, -1, __FUNCTION__, extra); +#endif + smpi_process_finalize(); +#ifdef HAVE_TRACING + TRACE_smpi_collective_out(rank, -1, __FUNCTION__); + TRACE_smpi_finalize(smpi_process_index()); +#endif + smpi_process_destroy(); + return MPI_SUCCESS; }