X-Git-Url: http://bilbo.iut-bm.univ-fcomte.fr/pub/gitweb/simgrid.git/blobdiff_plain/92959d5b26387e5194b1bb4553baf90da7701f41..401273041d310cd9077459a90fc18aa0361b034c:/src/smpi/smpi_replay.c diff --git a/src/smpi/smpi_replay.c b/src/smpi/smpi_replay.c index e41a4f29cd..9151825c05 100644 --- a/src/smpi/smpi_replay.c +++ b/src/smpi/smpi_replay.c @@ -170,7 +170,7 @@ static void action_compute(const char *const *action) double flops= parse_double(action[2]); #ifdef HAVE_TRACING int rank = smpi_comm_rank(MPI_COMM_WORLD); - instr_extra_data extra = xbt_new(s_instr_extra_data_t,1); + 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); @@ -199,7 +199,7 @@ 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); - instr_extra_data extra = xbt_new(s_instr_extra_data_t,1); + instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1); extra->type = TRACING_SEND; extra->send_size = size; extra->src = rank; @@ -232,7 +232,7 @@ 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); - instr_extra_data extra = xbt_new(s_instr_extra_data_t,1); + instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1); extra->type = TRACING_ISEND; extra->send_size = size; extra->src = rank; @@ -267,7 +267,7 @@ 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); - instr_extra_data extra = xbt_new(s_instr_extra_data_t,1); + instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1); extra->type = TRACING_RECV; extra->send_size = size; extra->src = src_traced; @@ -302,7 +302,7 @@ 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); - instr_extra_data extra = xbt_new(s_instr_extra_data_t,1); + instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1); extra->type = TRACING_IRECV; extra->send_size = size; extra->src = src_traced; @@ -343,7 +343,7 @@ 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; - instr_extra_data extra = xbt_new(s_instr_extra_data_t,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 @@ -402,7 +402,7 @@ static void action_waitall(const char *const *action){ } } int rank_traced = smpi_process_index(); - instr_extra_data extra = xbt_new(s_instr_extra_data_t,1); + 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); @@ -436,7 +436,7 @@ 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); - instr_extra_data extra = xbt_new(s_instr_extra_data_t,1); + 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 @@ -471,7 +471,7 @@ static void action_bcast(const char *const *action) int rank = smpi_comm_rank(MPI_COMM_WORLD); int root_traced = smpi_group_index(smpi_comm_group(MPI_COMM_WORLD), root); - instr_extra_data extra = xbt_new(s_instr_extra_data_t,1); + instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1); extra->type = TRACING_BCAST; extra->send_size = size; extra->root = root_traced; @@ -506,7 +506,7 @@ 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), root); - instr_extra_data extra = xbt_new(s_instr_extra_data_t,1); + 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; @@ -534,7 +534,7 @@ 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); - instr_extra_data extra = xbt_new(s_instr_extra_data_t,1); + 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; @@ -572,7 +572,7 @@ static void action_allToAll(const char *const *action) { #ifdef HAVE_TRACING int rank = smpi_process_index(); - instr_extra_data extra = xbt_new(s_instr_extra_data_t,1); + 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; @@ -620,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(); @@ -629,7 +629,7 @@ static void action_gather(const char *const *action) { recv = calloc(recv_size*comm_size, smpi_datatype_size(MPI_CURRENT_TYPE2)); #ifdef HAVE_TRACING - instr_extra_data extra = xbt_new(s_instr_extra_data_t,1); + 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; @@ -653,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) { /* @@ -690,7 +763,7 @@ static void action_reducescatter(const char *const *action) { } #ifdef HAVE_TRACING - instr_extra_data extra = xbt_new(s_instr_extra_data_t,1); + 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)); @@ -713,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); } @@ -760,7 +834,7 @@ static void action_allgatherv(const char *const *action) { #ifdef HAVE_TRACING int rank = smpi_process_index(); - instr_extra_data extra = xbt_new(s_instr_extra_data_t,1); + 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)); @@ -835,11 +909,8 @@ static void action_allToAllv(const char *const *action) { #ifdef HAVE_TRACING int rank = smpi_process_index(); - int count=0; - for(i=0;itype = TRACING_ALLTOALLV; - extra->send_size = count; extra->recvcounts= xbt_malloc(comm_size*sizeof(int)); extra->sendcounts= xbt_malloc(comm_size*sizeof(int)); extra->num_processes = comm_size; @@ -847,6 +918,7 @@ static void action_allToAllv(const char *const *action) { 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); @@ -877,7 +949,7 @@ void smpi_replay_init(int *argc, char***argv){ int rank = smpi_process_index(); TRACE_smpi_init(rank); TRACE_smpi_computing_init(rank); - instr_extra_data extra = xbt_new(s_instr_extra_data_t,1); + 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__); @@ -903,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); @@ -930,7 +1003,7 @@ int smpi_replay_finalize(){ smpi_mpi_barrier(MPI_COMM_WORLD); #ifdef HAVE_TRACING int rank = smpi_process_index(); - instr_extra_data extra = xbt_new(s_instr_extra_data_t,1); + 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