-/* Copyright (c) 2009-2013. The SimGrid Team.
+/* Copyright (c) 2009-2014. The SimGrid Team.
* All rights reserved. */
/* This program is free software; you can redistribute it and/or modify it
* under the terms of the license (GNU LGPL) which comes with this package. */
-
#include "private.h"
#include <stdio.h>
#include <xbt.h>
MPI_Datatype MPI_DEFAULT_TYPE;
MPI_Datatype MPI_CURRENT_TYPE;
+static int sendbuffer_size=0;
+char* sendbuffer=NULL;
+static int recvbuffer_size=0;
+char* recvbuffer=NULL;
+
static void log_timed_action (const char *const *action, double clock){
if (XBT_LOG_ISENABLED(smpi_replay, xbt_log_priority_verbose)){
char *name = xbt_str_join_array(action, " ");
}
}
-typedef struct {
- xbt_dynar_t irecvs; /* of MPI_Request */
-} s_smpi_replay_globals_t, *smpi_replay_globals_t;
+//allocate a single buffer for all sends, growing it if needed
+void* smpi_get_tmp_sendbuffer(int size){
+ if (!_xbt_replay_is_active())
+ return xbt_malloc(size);
+ if (sendbuffer_size<size){
+ sendbuffer=xbt_realloc(sendbuffer,size);
+ sendbuffer_size=size;
+ }
+ return sendbuffer;
+}
+//allocate a single buffer for all recv
+void* smpi_get_tmp_recvbuffer(int size){
+ if (!_xbt_replay_is_active())
+ return xbt_malloc(size);
+ if (recvbuffer_size<size){
+ recvbuffer=xbt_realloc(recvbuffer,size);
+ recvbuffer_size=size;
+ }
+ return sendbuffer;
+}
+void smpi_free_tmp_buffer(void* buf){
+ if (!_xbt_replay_is_active())
+ xbt_free(buf);
+}
/* Helper function */
static double parse_double(const char *string)
{
int i;
XBT_DEBUG("Initialize the counters");
- smpi_replay_globals_t globals = xbt_new(s_smpi_replay_globals_t, 1);
- globals->irecvs = xbt_dynar_new(sizeof(MPI_Request),NULL);
if(action[2]) MPI_DEFAULT_TYPE= MPI_DOUBLE; // default MPE dataype
else MPI_DEFAULT_TYPE= MPI_BYTE; // default TAU datatype
- smpi_process_set_user_data((void*) globals);
-
/* start a simulated timer */
smpi_process_simulated_start();
/*initialize the number of active processes */
reqq=xbt_new0(xbt_dynar_t,active_processes);
for(i=0;i<active_processes;i++){
- reqq[i]=xbt_dynar_new(sizeof(MPI_Request),NULL);
+ reqq[i]=xbt_dynar_new(sizeof(MPI_Request),&xbt_free_ref);
}
}
}
static void action_finalize(const char *const *action)
{
- smpi_replay_globals_t globals =
- (smpi_replay_globals_t) smpi_process_get_user_data();
- if (globals){
- XBT_DEBUG("There are %lu irecvs in the dynar",
- xbt_dynar_length(globals->irecvs));
- xbt_dynar_free_container(&(globals->irecvs));
- }
- free(globals);
}
static void action_comm_size(const char *const *action)
double clock = smpi_process_simulated_elapsed();
double flops= parse_double(action[2]);
#ifdef HAVE_TRACING
- int rank = smpi_comm_rank(MPI_COMM_WORLD);
+ int rank = smpi_process_index();
instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
extra->type=TRACING_COMPUTING;
extra->comp_size=flops;
}
#ifdef HAVE_TRACING
- int rank = smpi_comm_rank(MPI_COMM_WORLD);
+ int rank = smpi_process_index();
int dst_traced = smpi_group_rank(smpi_comm_group(MPI_COMM_WORLD), to);
instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
else MPI_CURRENT_TYPE= MPI_DEFAULT_TYPE;
#ifdef HAVE_TRACING
- int rank = smpi_comm_rank(MPI_COMM_WORLD);
+ int rank = smpi_process_index();
int dst_traced = smpi_group_rank(smpi_comm_group(MPI_COMM_WORLD), to);
instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
extra->type = TRACING_ISEND;
request->send = 1;
#endif
- xbt_dynar_push(reqq[smpi_comm_rank(MPI_COMM_WORLD)],&request);
+ xbt_dynar_push(reqq[smpi_process_index()],&request);
log_timed_action (action, clock);
}
else MPI_CURRENT_TYPE= MPI_DEFAULT_TYPE;
#ifdef HAVE_TRACING
- int rank = smpi_comm_rank(MPI_COMM_WORLD);
+ int rank = smpi_process_index();
int src_traced = smpi_group_rank(smpi_comm_group(MPI_COMM_WORLD), from);
instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__, extra);
#endif
+ //unknow size from the receiver pov
+ if(size==-1){
+ smpi_mpi_probe(from, 0, MPI_COMM_WORLD, &status);
+ size=status.count;
+ }
+
smpi_mpi_recv(NULL, size, MPI_CURRENT_TYPE, from, 0, MPI_COMM_WORLD, &status);
#ifdef HAVE_TRACING
double clock = smpi_process_simulated_elapsed();
MPI_Request request;
- smpi_replay_globals_t globals =
- (smpi_replay_globals_t) smpi_process_get_user_data();
-
if(action[4]) MPI_CURRENT_TYPE=decode_datatype(action[4]);
else MPI_CURRENT_TYPE= MPI_DEFAULT_TYPE;
#ifdef HAVE_TRACING
- int rank = smpi_comm_rank(MPI_COMM_WORLD);
+ int rank = smpi_process_index();
int src_traced = smpi_group_rank(smpi_comm_group(MPI_COMM_WORLD), from);
instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
extra->type = TRACING_IRECV;
extra->datatype1 = encode_datatype(MPI_CURRENT_TYPE);
TRACE_smpi_ptp_in(rank, src_traced, rank, __FUNCTION__, extra);
#endif
+ MPI_Status status;
+ //unknow size from the receiver pov
+ if(size==-1){
+ smpi_mpi_probe(from, 0, MPI_COMM_WORLD, &status);
+ size=status.count;
+ }
request = smpi_mpi_irecv(NULL, size, MPI_CURRENT_TYPE, from, 0, MPI_COMM_WORLD);
TRACE_smpi_ptp_out(rank, src_traced, rank, __FUNCTION__);
request->recv = 1;
#endif
- xbt_dynar_push(globals->irecvs,&request);
- xbt_dynar_push(reqq[smpi_comm_rank(MPI_COMM_WORLD)],&request);
+ xbt_dynar_push(reqq[smpi_process_index()],&request);
+
+ log_timed_action (action, clock);
+}
+
+static void action_test(const char *const *action){
+ double clock = smpi_process_simulated_elapsed();
+ MPI_Request request;
+ MPI_Status status;
+ int flag = TRUE;
+
+ request = xbt_dynar_pop_as(reqq[smpi_process_index()],MPI_Request);
+ xbt_assert(request != NULL, "found null request in reqq");
+
+#ifdef HAVE_TRACING
+ int rank = smpi_process_index();
+ instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
+ extra->type=TRACING_TEST;
+ TRACE_smpi_testing_in(rank, extra);
+#endif
+ flag = smpi_mpi_test(&request, &status);
+ XBT_DEBUG("MPI_Test result: %d", flag);
+ /* push back request in dynar to be caught by a subsequent wait. if the test
+ * did succeed, the request is now NULL.
+ */
+ xbt_dynar_push_as(reqq[smpi_process_index()],MPI_Request, request);
+
+#ifdef HAVE_TRACING
+ TRACE_smpi_testing_out(rank);
+#endif
log_timed_action (action, clock);
}
double clock = smpi_process_simulated_elapsed();
MPI_Request request;
MPI_Status status;
- smpi_replay_globals_t globals =
- (smpi_replay_globals_t) smpi_process_get_user_data();
- xbt_assert(xbt_dynar_length(globals->irecvs),
- "action wait not preceded by any irecv: %s",
+ xbt_assert(xbt_dynar_length(reqq[smpi_process_index()]),
+ "action wait not preceded by any irecv or isend: %s",
xbt_str_join_array(action," "));
- request = xbt_dynar_pop_as(globals->irecvs,MPI_Request);
- xbt_assert(request != NULL, "found null request in globals->irecv");
+ request = xbt_dynar_pop_as(reqq[smpi_process_index()],MPI_Request);
+
+ if (!request){
+ /* Assuming that the trace is well formed, this mean the comm might have
+ * been caught by a MPI_test. Then just return.
+ */
+ return;
+ }
+
#ifdef HAVE_TRACING
int rank = request->comm != MPI_COMM_NULL
? smpi_comm_rank(request->comm)
int count_requests=0;
unsigned int i=0;
- count_requests=xbt_dynar_length(reqq[smpi_comm_rank(MPI_COMM_WORLD)]);
+ count_requests=xbt_dynar_length(reqq[smpi_process_index()]);
if (count_requests>0) {
MPI_Request requests[count_requests];
/* The reqq is an array of dynars. Its index corresponds to the rank.
Thus each rank saves its own requests to the array request. */
- xbt_dynar_foreach(reqq[smpi_comm_rank(MPI_COMM_WORLD)],i,requests[i]);
+ xbt_dynar_foreach(reqq[smpi_process_index()],i,requests[i]);
#ifdef HAVE_TRACING
//save information from requests
xbt_dynar_free(&recvs);
#endif
- xbt_dynar_reset(reqq[smpi_comm_rank(MPI_COMM_WORLD)]);
+ xbt_dynar_free_container(&(reqq[smpi_process_index()]));
}
log_timed_action (action, clock);
}
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);
+ int rank = smpi_process_index();
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);
+ mpi_coll_barrier_fun(MPI_COMM_WORLD);
#ifdef HAVE_TRACING
TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
#endif
}
#ifdef HAVE_TRACING
- int rank = smpi_comm_rank(MPI_COMM_WORLD);
+ int rank = smpi_process_index();
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);
}
#ifdef HAVE_TRACING
- int rank = smpi_comm_rank(MPI_COMM_WORLD);
+ int rank = smpi_process_index();
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;
double clock = smpi_process_simulated_elapsed();
#ifdef HAVE_TRACING
- int rank = smpi_comm_rank(MPI_COMM_WORLD);
+ int rank = smpi_process_index();
instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
extra->type = TRACING_ALLREDUCE;
extra->send_size = comm_size;
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);
+ mpi_coll_allreduce_fun(NULL, NULL, comm_size, MPI_CURRENT_TYPE, MPI_OP_NULL, MPI_COMM_WORLD);
smpi_execute_flops(comp_size);
- mpi_coll_bcast_fun(NULL, comm_size, MPI_CURRENT_TYPE, 0, MPI_COMM_WORLD);
#ifdef HAVE_TRACING
TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
#endif
MPI_CURRENT_TYPE=MPI_DEFAULT_TYPE;
MPI_CURRENT_TYPE2=MPI_DEFAULT_TYPE;
}
- void *send = calloc(send_size*comm_size, smpi_datatype_size(MPI_CURRENT_TYPE));
- void *recv = calloc(recv_size*comm_size, smpi_datatype_size(MPI_CURRENT_TYPE2));
+ void *send = smpi_get_tmp_sendbuffer(send_size*comm_size* smpi_datatype_size(MPI_CURRENT_TYPE));
+ void *recv = smpi_get_tmp_recvbuffer(recv_size*comm_size* smpi_datatype_size(MPI_CURRENT_TYPE2));
#ifdef HAVE_TRACING
int rank = smpi_process_index();
#endif
log_timed_action (action, clock);
- xbt_free(send);
- xbt_free(recv);
+
}
MPI_CURRENT_TYPE=MPI_DEFAULT_TYPE;
MPI_CURRENT_TYPE2=MPI_DEFAULT_TYPE;
}
- void *send = calloc(send_size, smpi_datatype_size(MPI_CURRENT_TYPE));
+ void *send = smpi_get_tmp_sendbuffer(send_size* smpi_datatype_size(MPI_CURRENT_TYPE));
void *recv = NULL;
int root=atoi(action[4]);
int rank = smpi_process_index();
if(rank==root)
- recv = calloc(recv_size*comm_size, smpi_datatype_size(MPI_CURRENT_TYPE2));
+ recv = smpi_get_tmp_recvbuffer(recv_size*comm_size* smpi_datatype_size(MPI_CURRENT_TYPE2));
#ifdef HAVE_TRACING
instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
TRACE_smpi_collective_in(rank, root, __FUNCTION__, extra);
#endif
-smpi_mpi_gather(send, send_size, MPI_CURRENT_TYPE,
+ mpi_coll_gather_fun(send, send_size, MPI_CURRENT_TYPE,
recv, recv_size, MPI_CURRENT_TYPE2,
root, MPI_COMM_WORLD);
#endif
log_timed_action (action, clock);
- xbt_free(send);
- xbt_free(recv);
+
}
MPI_CURRENT_TYPE=MPI_DEFAULT_TYPE;
MPI_CURRENT_TYPE2=MPI_DEFAULT_TYPE;
}
- void *send = calloc(send_size, smpi_datatype_size(MPI_CURRENT_TYPE));
+ void *send = smpi_get_tmp_sendbuffer(send_size* smpi_datatype_size(MPI_CURRENT_TYPE));
void *recv = NULL;
for(i=0;i<comm_size;i++) {
recvcounts[i] = atoi(action[i+3]);
int rank = smpi_process_index();
if(rank==root)
- recv = calloc(recv_sum, smpi_datatype_size(MPI_CURRENT_TYPE2));
+ recv = smpi_get_tmp_recvbuffer(recv_sum* smpi_datatype_size(MPI_CURRENT_TYPE2));
#ifdef HAVE_TRACING
instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);
log_timed_action (action, clock);
xbt_free(recvcounts);
- xbt_free(send);
- xbt_free(recv);
xbt_free(disps);
}
int comp_size = parse_double(action[2+comm_size]);
int *recvcounts = xbt_new0(int, comm_size);
int *disps = xbt_new0(int, comm_size);
- int i=0,recv_sum=0;
- int root=0;
+ int i=0;
int rank = smpi_process_index();
if(action[3+comm_size])
for(i=0;i<comm_size;i++) {
recvcounts[i] = atoi(action[i+2]);
- recv_sum=recv_sum+recvcounts[i];
disps[i] = 0;
}
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);
- smpi_mpi_scatterv(NULL, recvcounts, disps, MPI_CURRENT_TYPE, NULL,
- recvcounts[rank], MPI_CURRENT_TYPE, 0, MPI_COMM_WORLD);
+ mpi_coll_reduce_scatter_fun(NULL, NULL, recvcounts, MPI_CURRENT_TYPE, MPI_OP_NULL,
+ MPI_COMM_WORLD);
smpi_execute_flops(comp_size);
MPI_CURRENT_TYPE = MPI_DEFAULT_TYPE;
MPI_CURRENT_TYPE2 = MPI_DEFAULT_TYPE;
}
- void *sendbuf = calloc(sendcount, smpi_datatype_size(MPI_CURRENT_TYPE));
+ void *sendbuf = smpi_get_tmp_sendbuffer(sendcount* smpi_datatype_size(MPI_CURRENT_TYPE));
for(i=0;i<comm_size;i++) {
recvcounts[i] = atoi(action[i+3]);
recv_sum=recv_sum+recvcounts[i];
}
- void *recvbuf = calloc(recv_sum, smpi_datatype_size(MPI_CURRENT_TYPE2));
+ void *recvbuf = smpi_get_tmp_recvbuffer(recv_sum* smpi_datatype_size(MPI_CURRENT_TYPE2));
#ifdef HAVE_TRACING
int rank = smpi_process_index();
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);
+ mpi_coll_allgatherv_fun(sendbuf, sendcount, MPI_CURRENT_TYPE, recvbuf, recvcounts, disps, MPI_CURRENT_TYPE2, MPI_COMM_WORLD);
#ifdef HAVE_TRACING
TRACE_smpi_collective_out(rank, -1, __FUNCTION__);
#endif
log_timed_action (action, clock);
- xbt_free(sendbuf);
- xbt_free(recvbuf);
xbt_free(recvcounts);
xbt_free(disps);
}
MPI_CURRENT_TYPE2=MPI_DEFAULT_TYPE;
}
- void *sendbuf = calloc(send_buf_size, smpi_datatype_size(MPI_CURRENT_TYPE));
- void *recvbuf = calloc(recv_buf_size, smpi_datatype_size(MPI_CURRENT_TYPE2));
+ void *sendbuf = smpi_get_tmp_sendbuffer(send_buf_size* smpi_datatype_size(MPI_CURRENT_TYPE));
+ void *recvbuf = smpi_get_tmp_recvbuffer(recv_buf_size* smpi_datatype_size(MPI_CURRENT_TYPE2));
for(i=0;i<comm_size;i++) {
sendcounts[i] = atoi(action[i+3]);
TRACE_smpi_collective_in(rank, -1, __FUNCTION__,extra);
#endif
- mpi_coll_alltoallv_fun(sendbuf, sendcounts, senddisps, MPI_CURRENT_TYPE,
+ mpi_coll_alltoallv_fun(sendbuf, sendcounts, senddisps, MPI_CURRENT_TYPE,
recvbuf, recvcounts, recvdisps, MPI_CURRENT_TYPE,
MPI_COMM_WORLD);
#ifdef HAVE_TRACING
#endif
log_timed_action (action, clock);
- xbt_free(sendbuf);
- xbt_free(recvbuf);
xbt_free(sendcounts);
xbt_free(recvcounts);
xbt_free(senddisps);
xbt_replay_action_register("Isend", action_Isend);
xbt_replay_action_register("recv", action_recv);
xbt_replay_action_register("Irecv", action_Irecv);
+ xbt_replay_action_register("test", action_test);
xbt_replay_action_register("wait", action_wait);
xbt_replay_action_register("waitAll", action_waitall);
xbt_replay_action_register("barrier", action_barrier);
xbt_replay_action_register("reduceScatter", action_reducescatter);
xbt_replay_action_register("compute", action_compute);
}
-
+
+ //if we have a delayed start, sleep here.
+ if(*argc>2){
+ char *endptr;
+ double value = strtod((*argv)[2], &endptr);
+ if (*endptr != '\0')
+ THROWF(unknown_error, 0, "%s is not a double", (*argv)[2]);
+ XBT_VERB("Delayed start for instance - Sleeping for %f flops ",value );
+ smpi_execute_flops(value);
+ }
xbt_replay_action_runner(*argc, *argv);
}
int smpi_replay_finalize(){
double sim_time= 1.;
/* One active process will stop. Decrease the counter*/
- active_processes--;
XBT_DEBUG("There are %lu elements in reqq[*]",
- xbt_dynar_length(reqq[smpi_comm_rank(MPI_COMM_WORLD)]));
- xbt_dynar_free(&reqq[smpi_comm_rank(MPI_COMM_WORLD)]);
+ xbt_dynar_length(reqq[smpi_process_index()]));
+ if (!xbt_dynar_is_empty(reqq[smpi_process_index()])){
+ int count_requests=xbt_dynar_length(reqq[smpi_process_index()]);
+ MPI_Request requests[count_requests];
+ MPI_Status status[count_requests];
+ unsigned int i;
+
+ xbt_dynar_foreach(reqq[smpi_process_index()],i,requests[i]);
+ smpi_mpi_waitall(count_requests, requests, status);
+ active_processes--;
+ } else {
+ active_processes--;
+ }
+
+ xbt_dynar_free_container(&(reqq[smpi_process_index()]));
+
if(!active_processes){
/* Last process alive speaking */
/* end the simulated timer */
sim_time = smpi_process_simulated_elapsed();
- XBT_INFO("Simulation time %g", sim_time);
+ XBT_INFO("Simulation time %f", sim_time);
_xbt_replay_action_exit();
+ xbt_free(sendbuffer);
+ xbt_free(recvbuffer);
xbt_free(reqq);
reqq = NULL;
}
- smpi_mpi_barrier(MPI_COMM_WORLD);
+ mpi_coll_barrier_fun(MPI_COMM_WORLD);
#ifdef HAVE_TRACING
int rank = smpi_process_index();
instr_extra_data extra = xbt_new0(s_instr_extra_data_t,1);