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
1 /* Copyright (c) 2007-2022. The SimGrid Team. All rights reserved.          */
2
3 /* This program is free software; you can redistribute it and/or modify it
4  * under the terms of the license (GNU LGPL) which comes with this package. */
5
6 #include "private.hpp"
7
8 #include "smpi_file.hpp"
9 #include "smpi_datatype.hpp"
10
11 XBT_LOG_EXTERNAL_DEFAULT_CATEGORY(smpi_pmpi);
12 #define CHECK_RDONLY(fh)                                                                                               \
13   if ((fh)->flags() & MPI_MODE_RDONLY)                                                                                 \
14     return MPI_ERR_AMODE;
15 #define CHECK_WRONLY(fh)                                                                                               \
16   if ((fh)->flags() & MPI_MODE_WRONLY)                                                                                 \
17     return MPI_ERR_AMODE;
18 #define PASS_ZEROCOUNT(count)                                                                                          \
19   if ((count) == 0) {                                                                                                  \
20     status->count = 0;                                                                                                 \
21     return MPI_SUCCESS;                                                                                                \
22   }
23
24 #define CHECK_FILE_INPUTS                                                                                              \
25   CHECK_FILE(1, fh)                                                                                                    \
26   CHECK_COUNT(3, count)                                                                                                \
27   CHECK_TYPE(4, datatype)                                                                                              \
28   CHECK_BUFFER(2, buf, count, datatype)                                                                                          \
29
30 #define CHECK_FILE_INPUT_OFFSET                                                                                        \
31   CHECK_FILE(1, fh)                                                                                                    \
32   CHECK_OFFSET(3, offset)                                                                                                 \
33   CHECK_COUNT(4, count)                                                                                                \
34   CHECK_TYPE(5, datatype)                                                                                              \
35   CHECK_BUFFER(2, buf, count, datatype)                                                                                          \
36
37 extern MPI_Errhandler SMPI_default_File_Errhandler;
38
39 int PMPI_File_open(MPI_Comm comm, const char *filename, int amode, MPI_Info info, MPI_File *fh){
40   CHECK_COMM(1)
41   CHECK_COLLECTIVE(comm, "MPI_File_open")
42   CHECK_NULL(2, MPI_ERR_FILE, filename)
43   if (amode < 0)
44     return MPI_ERR_AMODE;
45   const SmpiBenchGuard suspend_bench;
46   *fh =  new simgrid::smpi::File(comm, filename, amode, info);
47   if ((*fh)->size() != 0 && (amode & MPI_MODE_EXCL)){
48     delete fh;
49     return MPI_ERR_AMODE;
50   }
51   if(amode & MPI_MODE_APPEND)
52     (*fh)->seek(0,MPI_SEEK_END);
53   return MPI_SUCCESS;
54 }
55
56 int PMPI_File_close(MPI_File *fh){
57   CHECK_NULL(2, MPI_ERR_ARG, fh)
58   CHECK_COLLECTIVE((*fh)->comm(), __func__)
59   const SmpiBenchGuard suspend_bench;
60   int ret = simgrid::smpi::File::close(fh);
61   *fh = MPI_FILE_NULL;
62   return ret;
63 }
64
65 int PMPI_File_seek(MPI_File fh, MPI_Offset offset, int whence){
66   CHECK_FILE(1, fh)
67   const SmpiBenchGuard suspend_bench;
68   int ret = fh->seek(offset,whence);
69   return ret;
70 }
71
72 int PMPI_File_seek_shared(MPI_File fh, MPI_Offset offset, int whence){
73   CHECK_FILE(1, fh)
74   CHECK_COLLECTIVE(fh->comm(), __func__)
75   const SmpiBenchGuard suspend_bench;
76   int ret = fh->seek_shared(offset,whence);
77   return ret;
78 }
79
80 int PMPI_File_get_position(MPI_File fh, MPI_Offset* offset){
81   CHECK_FILE(1, fh)
82   CHECK_NULL(2, MPI_ERR_DISP, offset)
83   const SmpiBenchGuard suspend_bench;
84   int ret = fh->get_position(offset);
85   return ret;
86 }
87
88 int PMPI_File_get_position_shared(MPI_File fh, MPI_Offset* offset){
89   CHECK_FILE(1, fh)
90   CHECK_NULL(2, MPI_ERR_DISP, offset)
91   const SmpiBenchGuard suspend_bench;
92   int ret = fh->get_position_shared(offset);
93   return ret;
94 }
95
96 int PMPI_File_read(MPI_File fh, void *buf, int count,MPI_Datatype datatype, MPI_Status *status){
97   CHECK_FILE_INPUTS
98   CHECK_WRONLY(fh)
99   PASS_ZEROCOUNT(count)
100   const SmpiBenchGuard suspend_bench;
101   aid_t rank_traced = simgrid::s4u::this_actor::get_pid();
102   TRACE_smpi_comm_in(rank_traced, __func__, new simgrid::instr::CpuTIData("IO - read", count * datatype->size()));
103   int ret = simgrid::smpi::File::read(fh, buf, count, datatype, status);
104   TRACE_smpi_comm_out(rank_traced);
105   return ret;
106 }
107
108 int PMPI_File_read_shared(MPI_File fh, void *buf, int count,MPI_Datatype datatype, MPI_Status *status){
109   CHECK_FILE_INPUTS
110   CHECK_WRONLY(fh)
111   PASS_ZEROCOUNT(count)
112   const SmpiBenchGuard suspend_bench;
113   aid_t rank_traced = simgrid::s4u::this_actor::get_pid();
114   TRACE_smpi_comm_in(rank_traced, __func__,
115                      new simgrid::instr::CpuTIData("IO - read_shared", count * datatype->size()));
116   int ret = simgrid::smpi::File::read_shared(fh, buf, count, datatype, status);
117   TRACE_smpi_comm_out(rank_traced);
118   return ret;
119 }
120
121 int PMPI_File_write(MPI_File fh, const void *buf, int count,MPI_Datatype datatype, MPI_Status *status){
122   CHECK_FILE_INPUTS
123   CHECK_RDONLY(fh)
124   PASS_ZEROCOUNT(count)
125   const SmpiBenchGuard suspend_bench;
126   aid_t rank_traced = simgrid::s4u::this_actor::get_pid();
127   TRACE_smpi_comm_in(rank_traced, __func__, new simgrid::instr::CpuTIData("IO - write", count * datatype->size()));
128   int ret = simgrid::smpi::File::write(fh, const_cast<void*>(buf), count, datatype, status);
129   TRACE_smpi_comm_out(rank_traced);
130   return ret;
131 }
132
133 int PMPI_File_write_shared(MPI_File fh, const void *buf, int count,MPI_Datatype datatype, MPI_Status *status){
134   CHECK_FILE_INPUTS
135   CHECK_RDONLY(fh)
136   PASS_ZEROCOUNT(count)
137   const SmpiBenchGuard suspend_bench;
138   aid_t rank_traced = simgrid::s4u::this_actor::get_pid();
139   TRACE_smpi_comm_in(rank_traced, __func__,
140                      new simgrid::instr::CpuTIData("IO - write_shared", count * datatype->size()));
141   int ret = simgrid::smpi::File::write_shared(fh, buf, count, datatype, status);
142   TRACE_smpi_comm_out(rank_traced);
143   return ret;
144 }
145
146 int PMPI_File_read_all(MPI_File fh, void *buf, int count,MPI_Datatype datatype, MPI_Status *status){
147   CHECK_FILE_INPUTS
148   CHECK_WRONLY(fh)
149   CHECK_COLLECTIVE(fh->comm(), __func__)
150   const SmpiBenchGuard suspend_bench;
151   aid_t rank_traced = simgrid::s4u::this_actor::get_pid();
152   TRACE_smpi_comm_in(rank_traced, __func__, new simgrid::instr::CpuTIData("IO - read_all", count * datatype->size()));
153   int ret = fh->op_all<simgrid::smpi::File::read>(buf, count, datatype, status);
154   TRACE_smpi_comm_out(rank_traced);
155   return ret;
156 }
157
158 int PMPI_File_read_ordered(MPI_File fh, void *buf, int count,MPI_Datatype datatype, MPI_Status *status){
159   CHECK_FILE_INPUTS
160   CHECK_WRONLY(fh)
161   CHECK_COLLECTIVE(fh->comm(), __func__)
162   const SmpiBenchGuard suspend_bench;
163   aid_t rank_traced = simgrid::s4u::this_actor::get_pid();
164   TRACE_smpi_comm_in(rank_traced, __func__,
165                      new simgrid::instr::CpuTIData("IO - read_ordered", count * datatype->size()));
166   int ret = simgrid::smpi::File::read_ordered(fh, buf, count, datatype, status);
167   TRACE_smpi_comm_out(rank_traced);
168   return ret;
169 }
170
171 int PMPI_File_write_all(MPI_File fh, const void *buf, int count,MPI_Datatype datatype, MPI_Status *status){
172   CHECK_FILE_INPUTS
173   CHECK_RDONLY(fh)
174   CHECK_COLLECTIVE(fh->comm(), __func__)
175   const SmpiBenchGuard suspend_bench;
176   aid_t rank_traced = simgrid::s4u::this_actor::get_pid();
177   TRACE_smpi_comm_in(rank_traced, __func__, new simgrid::instr::CpuTIData("IO - write_all", count * datatype->size()));
178   int ret = fh->op_all<simgrid::smpi::File::write>(const_cast<void*>(buf), count, datatype, status);
179   TRACE_smpi_comm_out(rank_traced);
180   return ret;
181 }
182
183 int PMPI_File_write_ordered(MPI_File fh, const void *buf, int count,MPI_Datatype datatype, MPI_Status *status){
184   CHECK_FILE_INPUTS
185   CHECK_RDONLY(fh)
186   CHECK_COLLECTIVE(fh->comm(), __func__)
187   const SmpiBenchGuard suspend_bench;
188   aid_t rank_traced = simgrid::s4u::this_actor::get_pid();
189   TRACE_smpi_comm_in(rank_traced, __func__,
190                      new simgrid::instr::CpuTIData("IO - write_ordered", count * datatype->size()));
191   int ret = simgrid::smpi::File::write_ordered(fh, buf, count, datatype, status);
192   TRACE_smpi_comm_out(rank_traced);
193   return ret;
194 }
195
196 int PMPI_File_read_at(MPI_File fh, MPI_Offset offset, void *buf, int count,MPI_Datatype datatype, MPI_Status *status){
197   CHECK_FILE_INPUTS
198   CHECK_WRONLY(fh)
199   PASS_ZEROCOUNT(count)
200   const SmpiBenchGuard suspend_bench;
201   aid_t rank_traced = simgrid::s4u::this_actor::get_pid();
202   TRACE_smpi_comm_in(rank_traced, __func__, new simgrid::instr::CpuTIData("IO - read", count * datatype->size()));
203   int ret = fh->seek(offset,MPI_SEEK_SET);
204   if (ret == MPI_SUCCESS)
205     ret = simgrid::smpi::File::read(fh, buf, count, datatype, status);
206   TRACE_smpi_comm_out(rank_traced);
207   return ret;
208 }
209
210 int PMPI_File_read_at_all(MPI_File fh, MPI_Offset offset, void *buf, int count,MPI_Datatype datatype, MPI_Status *status){
211   CHECK_FILE_INPUT_OFFSET
212   CHECK_WRONLY(fh)
213   CHECK_COLLECTIVE(fh->comm(), __func__)
214   const SmpiBenchGuard suspend_bench;
215   aid_t rank_traced = simgrid::s4u::this_actor::get_pid();
216   TRACE_smpi_comm_in(rank_traced, __func__,
217                      new simgrid::instr::CpuTIData("IO - read_at_all", count * datatype->size()));
218   int ret = fh->seek(offset,MPI_SEEK_SET);
219   if (ret == MPI_SUCCESS)
220     ret = fh->op_all<simgrid::smpi::File::read>(buf, count, datatype, status);
221   TRACE_smpi_comm_out(rank_traced);
222   return ret;
223 }
224
225 int PMPI_File_write_at(MPI_File fh, MPI_Offset offset, const void *buf, int count,MPI_Datatype datatype, MPI_Status *status){
226   CHECK_FILE_INPUT_OFFSET
227   CHECK_RDONLY(fh)
228   PASS_ZEROCOUNT(count)
229   const SmpiBenchGuard suspend_bench;
230   aid_t rank_traced = simgrid::s4u::this_actor::get_pid();
231   TRACE_smpi_comm_in(rank_traced, __func__, new simgrid::instr::CpuTIData("IO - write", count * datatype->size()));
232   int ret = fh->seek(offset,MPI_SEEK_SET);
233   if (ret == MPI_SUCCESS)
234     ret = simgrid::smpi::File::write(fh, const_cast<void*>(buf), count, datatype, status);
235   TRACE_smpi_comm_out(rank_traced);
236   return ret;
237 }
238
239 int PMPI_File_write_at_all(MPI_File fh, MPI_Offset offset, const void *buf, int count,MPI_Datatype datatype, MPI_Status *status){
240   CHECK_FILE_INPUT_OFFSET
241   CHECK_RDONLY(fh)
242   CHECK_COLLECTIVE(fh->comm(), __func__)
243   const SmpiBenchGuard suspend_bench;
244   aid_t rank_traced = simgrid::s4u::this_actor::get_pid();
245   TRACE_smpi_comm_in(rank_traced, __func__,
246                      new simgrid::instr::CpuTIData("IO - write_at_all", count * datatype->size()));
247   int ret = fh->seek(offset,MPI_SEEK_SET);
248   if (ret == MPI_SUCCESS)
249     ret = fh->op_all<simgrid::smpi::File::write>(const_cast<void*>(buf), count, datatype, status);
250   TRACE_smpi_comm_out(rank_traced);
251   return ret;
252 }
253
254 int PMPI_File_delete(const char *filename, MPI_Info info){
255   if (filename == nullptr)
256     return MPI_ERR_FILE;
257   const SmpiBenchGuard suspend_bench;
258   int ret = simgrid::smpi::File::del(filename, info);
259   return ret;
260 }
261
262 int PMPI_File_set_view(MPI_File fh, MPI_Offset disp, MPI_Datatype etype, MPI_Datatype filetype, const char *datarep, MPI_Info info){
263   CHECK_FILE(1, fh)
264   CHECK_COLLECTIVE(fh->comm(), __func__)
265   if(not ((fh->flags() & MPI_MODE_SEQUENTIAL) && (disp == MPI_DISPLACEMENT_CURRENT)))
266     CHECK_OFFSET(2, disp)
267   CHECK_TYPE(3, etype)
268   CHECK_TYPE(4, filetype)
269   const SmpiBenchGuard suspend_bench;
270   int ret = fh->set_view(disp, etype, filetype, datarep, info);
271   return ret;
272 }
273
274 int PMPI_File_get_view(MPI_File fh, MPI_Offset *disp, MPI_Datatype *etype, MPI_Datatype *filetype, char *datarep){
275   CHECK_FILE(1, fh)
276   CHECK_NULL(2, MPI_ERR_ARG, disp)
277   CHECK_NULL(3, MPI_ERR_ARG, etype)
278   CHECK_NULL(4, MPI_ERR_ARG, filetype)
279   const SmpiBenchGuard suspend_bench;
280   int ret = fh->get_view(disp, etype, filetype, datarep);
281   return ret;
282 }
283
284 int PMPI_File_get_info(MPI_File  fh, MPI_Info* info)
285 {
286   CHECK_FILE(1, fh)
287   *info = new simgrid::smpi::Info(fh->info());
288   return MPI_SUCCESS;
289 }
290
291 int PMPI_File_set_info(MPI_File  fh, MPI_Info info)
292 {
293   CHECK_FILE(1, fh)
294   fh->set_info(info);
295   return MPI_SUCCESS;
296 }
297
298 int PMPI_File_get_size(MPI_File  fh, MPI_Offset* size)
299 {
300   CHECK_FILE(1, fh)
301   *size = fh->size();
302   return MPI_SUCCESS;
303 }
304
305 int PMPI_File_set_size(MPI_File  fh, MPI_Offset size)
306 {
307   CHECK_FILE(1, fh)
308   CHECK_COLLECTIVE(fh->comm(), __func__)
309   fh->set_size(size);
310   return MPI_SUCCESS;
311 }
312
313 int PMPI_File_get_amode(MPI_File  fh, int* amode)
314 {
315   CHECK_FILE(1, fh)
316   *amode = fh->flags();
317   return MPI_SUCCESS;
318 }
319
320 int PMPI_File_get_group(MPI_File  fh, MPI_Group* group)
321 {
322   CHECK_FILE(1, fh)
323   *group = fh->comm()->group();
324   return MPI_SUCCESS;
325 }
326
327 int PMPI_File_sync(MPI_File  fh)
328 {
329   CHECK_FILE(1, fh)
330   fh->sync();
331   return MPI_SUCCESS;
332 }
333
334 int PMPI_File_create_errhandler(MPI_File_errhandler_function* function, MPI_Errhandler* errhandler){
335   *errhandler=new simgrid::smpi::Errhandler(function);
336   return MPI_SUCCESS;
337 }
338
339 int PMPI_File_get_errhandler(MPI_File file, MPI_Errhandler* errhandler){
340   if (errhandler==nullptr){
341     return MPI_ERR_ARG;
342   } else if (file == MPI_FILE_NULL) {
343     *errhandler = SMPI_default_File_Errhandler;
344     return MPI_SUCCESS;
345   }
346   *errhandler=file->errhandler();
347   return MPI_SUCCESS;
348 }
349
350 int PMPI_File_set_errhandler(MPI_File file, MPI_Errhandler errhandler){
351   if (errhandler==nullptr){
352     return MPI_ERR_ARG;
353   } else if (file == MPI_FILE_NULL) {
354     SMPI_default_File_Errhandler = errhandler;
355     return MPI_SUCCESS;
356   }
357   file->set_errhandler(errhandler);
358   return MPI_SUCCESS;
359 }
360
361 int PMPI_File_call_errhandler(MPI_File file,int errorcode){
362   if (file == nullptr) {
363     return MPI_ERR_WIN;
364   }
365   MPI_Errhandler err = file->errhandler();
366   err->call(file, errorcode);
367   simgrid::smpi::Errhandler::unref(err);
368   return MPI_SUCCESS;
369 }