Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
add check for collectives, using check_collectives_ordering utility.
[simgrid.git] / src / smpi / bindings / smpi_pmpi_file.cpp
index 70df57f..e87680f 100644 (file)
@@ -38,6 +38,7 @@ extern MPI_Errhandler SMPI_default_File_Errhandler;
 
 int PMPI_File_open(MPI_Comm comm, const char *filename, int amode, MPI_Info info, MPI_File *fh){
   CHECK_COMM(1)
+  CHECK_COLLECTIVE(comm, "MPI_File_open")
   CHECK_NULL(2, MPI_ERR_FILE, filename)
   if (amode < 0)
     return MPI_ERR_AMODE;
@@ -54,6 +55,7 @@ int PMPI_File_open(MPI_Comm comm, const char *filename, int amode, MPI_Info info
 
 int PMPI_File_close(MPI_File *fh){
   CHECK_NULL(2, MPI_ERR_ARG, fh)
+  CHECK_COLLECTIVE((*fh)->comm(), __func__)
   const SmpiBenchGuard suspend_bench;
   int ret = simgrid::smpi::File::close(fh);
   *fh = MPI_FILE_NULL;
@@ -69,6 +71,7 @@ int PMPI_File_seek(MPI_File fh, MPI_Offset offset, int whence){
 
 int PMPI_File_seek_shared(MPI_File fh, MPI_Offset offset, int whence){
   CHECK_FILE(1, fh)
+  CHECK_COLLECTIVE(fh->comm(), __func__)
   const SmpiBenchGuard suspend_bench;
   int ret = fh->seek_shared(offset,whence);
   return ret;
@@ -143,6 +146,7 @@ int PMPI_File_write_shared(MPI_File fh, const void *buf, int count,MPI_Datatype
 int PMPI_File_read_all(MPI_File fh, void *buf, int count,MPI_Datatype datatype, MPI_Status *status){
   CHECK_FILE_INPUTS
   CHECK_WRONLY(fh)
+  CHECK_COLLECTIVE(fh->comm(), __func__)
   const SmpiBenchGuard suspend_bench;
   aid_t rank_traced = simgrid::s4u::this_actor::get_pid();
   TRACE_smpi_comm_in(rank_traced, __func__, new simgrid::instr::CpuTIData("IO - read_all", count * datatype->size()));
@@ -154,6 +158,7 @@ int PMPI_File_read_all(MPI_File fh, void *buf, int count,MPI_Datatype datatype,
 int PMPI_File_read_ordered(MPI_File fh, void *buf, int count,MPI_Datatype datatype, MPI_Status *status){
   CHECK_FILE_INPUTS
   CHECK_WRONLY(fh)
+  CHECK_COLLECTIVE(fh->comm(), __func__)
   const SmpiBenchGuard suspend_bench;
   aid_t rank_traced = simgrid::s4u::this_actor::get_pid();
   TRACE_smpi_comm_in(rank_traced, __func__,
@@ -166,6 +171,7 @@ int PMPI_File_read_ordered(MPI_File fh, void *buf, int count,MPI_Datatype dataty
 int PMPI_File_write_all(MPI_File fh, const void *buf, int count,MPI_Datatype datatype, MPI_Status *status){
   CHECK_FILE_INPUTS
   CHECK_RDONLY(fh)
+  CHECK_COLLECTIVE(fh->comm(), __func__)
   const SmpiBenchGuard suspend_bench;
   aid_t rank_traced = simgrid::s4u::this_actor::get_pid();
   TRACE_smpi_comm_in(rank_traced, __func__, new simgrid::instr::CpuTIData("IO - write_all", count * datatype->size()));
@@ -177,6 +183,7 @@ int PMPI_File_write_all(MPI_File fh, const void *buf, int count,MPI_Datatype dat
 int PMPI_File_write_ordered(MPI_File fh, const void *buf, int count,MPI_Datatype datatype, MPI_Status *status){
   CHECK_FILE_INPUTS
   CHECK_RDONLY(fh)
+  CHECK_COLLECTIVE(fh->comm(), __func__)
   const SmpiBenchGuard suspend_bench;
   aid_t rank_traced = simgrid::s4u::this_actor::get_pid();
   TRACE_smpi_comm_in(rank_traced, __func__,
@@ -203,6 +210,7 @@ int PMPI_File_read_at(MPI_File fh, MPI_Offset offset, void *buf, int count,MPI_D
 int PMPI_File_read_at_all(MPI_File fh, MPI_Offset offset, void *buf, int count,MPI_Datatype datatype, MPI_Status *status){
   CHECK_FILE_INPUT_OFFSET
   CHECK_WRONLY(fh)
+  CHECK_COLLECTIVE(fh->comm(), __func__)
   const SmpiBenchGuard suspend_bench;
   aid_t rank_traced = simgrid::s4u::this_actor::get_pid();
   TRACE_smpi_comm_in(rank_traced, __func__,
@@ -231,6 +239,7 @@ int PMPI_File_write_at(MPI_File fh, MPI_Offset offset, const void *buf, int coun
 int PMPI_File_write_at_all(MPI_File fh, MPI_Offset offset, const void *buf, int count,MPI_Datatype datatype, MPI_Status *status){
   CHECK_FILE_INPUT_OFFSET
   CHECK_RDONLY(fh)
+  CHECK_COLLECTIVE(fh->comm(), __func__)
   const SmpiBenchGuard suspend_bench;
   aid_t rank_traced = simgrid::s4u::this_actor::get_pid();
   TRACE_smpi_comm_in(rank_traced, __func__,
@@ -252,6 +261,7 @@ int PMPI_File_delete(const char *filename, MPI_Info info){
 
 int PMPI_File_set_view(MPI_File fh, MPI_Offset disp, MPI_Datatype etype, MPI_Datatype filetype, const char *datarep, MPI_Info info){
   CHECK_FILE(1, fh)
+  CHECK_COLLECTIVE(fh->comm(), __func__)
   if(not ((fh->flags() & MPI_MODE_SEQUENTIAL) && (disp == MPI_DISPLACEMENT_CURRENT)))
     CHECK_OFFSET(2, disp)
   CHECK_TYPE(3, etype)
@@ -295,6 +305,7 @@ int PMPI_File_get_size(MPI_File  fh, MPI_Offset* size)
 int PMPI_File_set_size(MPI_File  fh, MPI_Offset size)
 {
   CHECK_FILE(1, fh)
+  CHECK_COLLECTIVE(fh->comm(), __func__)
   fh->set_size(size);
   return MPI_SUCCESS;
 }