1 /*============================================================================
2  * Unit test for fvm_file.c;
3  *============================================================================*/
4 
5 /*
6   This file is part of Code_Saturne, a general-purpose CFD tool.
7 
8   Copyright (C) 1998-2021 EDF S.A.
9 
10   This program is free software; you can redistribute it and/or modify it under
11   the terms of the GNU General Public License as published by the Free Software
12   Foundation; either version 2 of the License, or (at your option) any later
13   version.
14 
15   This program is distributed in the hope that it will be useful, but WITHOUT
16   ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
17   FOR A PARTICULAR PURPOSE.  See the GNU General Public License for more
18   details.
19 
20   You should have received a copy of the GNU General Public License along with
21   this program; if not, write to the Free Software Foundation, Inc., 51 Franklin
22   Street, Fifth Floor, Boston, MA 02110-1301, USA.
23 */
24 
25 /*----------------------------------------------------------------------------*/
26 
27 #include "cs_defs.h"
28 
29 #include <stdio.h>
30 #include <stdlib.h>
31 #include <string.h>
32 
33 #include <bft_error.h>
34 #include <bft_mem.h>
35 #include <bft_printf.h>
36 
37 #include "cs_file.h"
38 
39 /*---------------------------------------------------------------------------*/
40 
41 static void
_create_test_data(void)42 _create_test_data(void)
43 {
44   int i;
45   cs_file_t *f;
46 
47   char header[80];
48   char footer[80];
49   int iarray[30];
50   double farray[30];
51 
52 #if defined(HAVE_MPI)
53   int mpi_flag;
54   MPI_Comm comm = MPI_COMM_NULL;
55 #endif
56 
57   sprintf(header, "fvm test file");
58   for (i = strlen(header); i < 80; i++)
59     header[i] = '\0';
60 
61   sprintf(footer, "fvm test file end");
62   for (i = strlen(footer); i < 80; i++)
63     footer[i] = '\0';
64 
65   for (i = 0; i < 30; i++)
66     iarray[i] = i+1;
67   for (i = 0; i < 30; i++)
68     farray[i] = i+1;
69 
70 #if defined(HAVE_MPI)
71   MPI_Initialized(&mpi_flag);
72   if (mpi_flag != 0) {
73     comm = MPI_COMM_WORLD;
74   }
75   f = cs_file_open("file_test_data",
76                    CS_FILE_MODE_WRITE,
77                    CS_FILE_STDIO_SERIAL,
78                    MPI_INFO_NULL,
79                    comm,
80                    comm);
81 #else
82   f = cs_file_open("file_test_data",
83                    CS_FILE_MODE_WRITE,
84                    0);
85 #endif
86 
87   cs_file_set_big_endian(f);
88 
89   cs_file_write_global(f, header, 1, 80);
90   cs_file_write_global(f, iarray, sizeof(int), 30);
91   cs_file_write_global(f, farray, sizeof(double), 30);
92 
93   cs_file_write_global(f, footer, 1, 80);
94 
95   f = cs_file_free(f);
96 }
97 
98 /*---------------------------------------------------------------------------*/
99 
100 int
main(int argc,char * argv[])101 main (int argc, char *argv[])
102 {
103   char mem_trace_name[32];
104   char output_file_name[32];
105   char buf[80];
106   int ibuf[30];
107   double dbuf[30];
108   cs_gnum_t i;
109   int a_id, p_id;
110   int size = 1;
111   int rank = 0;
112   size_t retval = 0;
113   cs_gnum_t block_start, block_end;
114   cs_gnum_t block_start_2, block_end_2;
115   cs_file_off_t off1 = -1, off2 = -1;
116   cs_file_t *f = NULL;
117 
118 #if defined(HAVE_MPI_IO)
119   const int n_pos = 2;
120   const int n_access = 5;
121   const cs_file_access_t access[5] = {CS_FILE_STDIO_SERIAL,
122                                       CS_FILE_STDIO_PARALLEL,
123                                       CS_FILE_MPI_INDEPENDENT,
124                                       CS_FILE_MPI_NON_COLLECTIVE,
125                                       CS_FILE_MPI_COLLECTIVE};
126   const cs_file_mpi_positioning_t pos[2] = {CS_FILE_MPI_EXPLICIT_OFFSETS,
127                                              CS_FILE_MPI_INDIVIDUAL_POINTERS};
128 #else
129   const int n_pos = 1;
130   const int n_access = 1;
131   const int access[1] = {CS_FILE_STDIO_SERIAL};
132   const cs_file_mpi_positioning_t pos[1] = {CS_FILE_MPI_EXPLICIT_OFFSETS};
133 #endif
134 
135 #if defined(HAVE_MPI)
136 
137   MPI_Status status;
138   int sync = 1;
139 
140   /* Initialization */
141 
142   MPI_Init(&argc, &argv);
143 
144   cs_glob_mpi_comm = MPI_COMM_WORLD;
145 
146   MPI_Comm_rank(MPI_COMM_WORLD, &rank);
147   MPI_Comm_size(MPI_COMM_WORLD, &size);
148 
149   block_start = rank * (30./size) + 1;
150   block_end   = (rank + 1) * (30./size) + 1;
151   if (rank == size - 1)
152     block_end = 31;
153 
154   block_start_2 = rank * (15./size) + 1;
155   block_end_2   = (rank + 1) * (15./size) + 1;
156   if (rank == size - 1)
157     block_end_2 = 16;
158 
159 #else
160 
161   block_start = 1;
162   block_end = 31;
163 
164   block_start_2 = 1;
165   block_end_2 = 16;
166 
167 #endif /* (HAVE_MPI) */
168 
169   if (size > 1)
170     sprintf(mem_trace_name, "cs_file_test_mem.%d", rank);
171   else
172     strcpy(mem_trace_name, "cs_file_test_mem");
173   bft_mem_init(mem_trace_name);
174 
175   _create_test_data();
176 
177   /* Loop on tests */
178 
179   for (a_id = 0; a_id < n_access; a_id++) {
180 
181     for (p_id = 0; p_id < n_pos; p_id++) {
182 
183       if (access[a_id] >= CS_FILE_MPI_INDEPENDENT) {
184 
185         cs_file_set_mpi_io_positioning(pos[p_id]);
186 
187         if (rank == 0)
188           bft_printf("Running test: %d-%d\n"
189                      "-------------\n\n", a_id, p_id);
190 
191         sprintf(output_file_name, "output_data_%d_%d", a_id+1, p_id+1);
192 
193       }
194       else {
195 
196         if (rank == 0)
197           bft_printf("Running test: %d\n"
198                      "-------------\n\n", a_id);
199 
200         sprintf(output_file_name, "output_data_%d", a_id+1);
201       }
202 
203       /* Read and seek/set tests */
204       /*-------------------------*/
205 
206 #if defined(HAVE_MPI)
207 
208       f = cs_file_open("file_test_data",
209                        CS_FILE_MODE_READ,
210                        access[a_id],
211                        MPI_INFO_NULL,
212                        MPI_COMM_WORLD,
213                        MPI_COMM_WORLD);
214 
215 #else
216 
217       f = cs_file_open("file_test_data",
218                        CS_FILE_MODE_READ,
219                        access[a_id]);
220 
221 #endif /* (HAVE_MPI) */
222 
223       cs_file_set_big_endian(f);
224 
225       cs_file_dump(f);
226 
227       retval = cs_file_read_global(f, buf, 1, 80);
228 
229       bft_printf("rank %d, readbuf = %s (returned %d)\n\n",
230                  rank, buf, (int)retval);
231 
232       for (i = 0; i < 30; i++)
233         ibuf[i] = 0;
234 
235       retval = cs_file_read_block(f, ibuf, sizeof(int), 1,
236                                   block_start, block_end);
237 
238 #if defined(HAVE_MPI) /* Serialize dump */
239       if (rank > 0)
240         MPI_Recv(&sync, 1, MPI_INT, rank - 1, 9876, MPI_COMM_WORLD, &status);
241 #endif
242 
243       bft_printf("\nRead by rank %d (returned %d):\n\n", rank, (int)retval);
244       for (i = block_start; i < block_end; i++)
245         bft_printf("  ival[%d] = %d\n", (int)i, (int)(ibuf[i-block_start]));
246 
247 #if defined(HAVE_MPI)
248       if (rank < size - 1)
249         MPI_Send(&sync, 1, MPI_INT, rank + 1, 9876, MPI_COMM_WORLD);
250 #endif
251 
252       off1 = cs_file_tell(f);
253 
254       retval = cs_file_read_block(f, dbuf, sizeof(double), 2,
255                                   block_start_2, block_end_2);
256 
257       off2 = cs_file_tell(f);
258 
259 #if defined(HAVE_MPI) /* Serialize dump */
260       if (rank > 0)
261         MPI_Recv(&sync, 1, MPI_INT, rank - 1, 9876, MPI_COMM_WORLD, &status);
262 #endif
263 
264       bft_printf("\nRead by rank %d (returned %d):\n\n", rank, (int)retval);
265       for (i = block_start_2; i < block_end_2; i++) {
266         bft_printf("  dval[%d] = %f\n", (int)(i*2 - 1),
267                    (dbuf[(i-block_start_2)*2]));
268         bft_printf("  dval[%d] = %f\n", (int)(i*2),
269                    (dbuf[(i-block_start_2)*2+1]));
270       }
271 
272 #if defined(HAVE_MPI)
273       if (rank < size - 1)
274         MPI_Send(&sync, 1, MPI_INT, rank + 1, 9876, MPI_COMM_WORLD);
275 #endif
276 
277 #if defined(HAVE_MPI)
278       MPI_Barrier(MPI_COMM_WORLD);
279 #endif
280       bft_printf("barrier passed by rank %d\n", rank);
281 
282       bft_printf("\nOffsets saved by rank %d: %ld, %ld\n\n",
283                  rank, (long)off1, (long)off2);
284 
285       retval = cs_file_read_global(f, buf, 1, 80);
286       bft_printf("rank %d, buf = %s (returned %d)\n", rank, buf, (int)retval);
287 
288       /* Test seek by re-reading at saved offset */
289 
290       cs_file_seek(f, off1, CS_FILE_SEEK_SET);
291 
292       memset(dbuf, 0, (block_end_2 - block_start_2)*sizeof(double));
293       retval = cs_file_read_block(f, dbuf, sizeof(double), 1,
294                                   block_start_2, block_end_2);
295 
296 #if defined(HAVE_MPI) /* Serialize dump */
297       if (rank > 0)
298         MPI_Recv(&sync, 1, MPI_INT, rank - 1, 9876, MPI_COMM_WORLD, &status);
299 #endif
300 
301       bft_printf("\nRe-read by rank %d (returned %d):\n\n", rank, (int)retval);
302       for (i = block_start; i < block_end; i++)
303         bft_printf("  dval[%d] = %f\n", (int)i, (dbuf[i-block_start]));
304 
305 #if defined(HAVE_MPI)
306       if (rank < size - 1)
307         MPI_Send(&sync, 1, MPI_INT, rank + 1, 9876, MPI_COMM_WORLD);
308 #endif
309 
310       cs_file_seek(f, off2, CS_FILE_SEEK_SET);
311 
312       retval = cs_file_read_global(f, buf, 1, 80);
313       bft_printf("rank %d, re-read buf = %s (returned %d)\n",
314                  rank, buf, (int)retval);
315 
316       f = cs_file_free(f);
317 
318       /* Write tests */
319       /*-------------*/
320 
321 #if defined(HAVE_MPI)
322 
323       f = cs_file_open(output_file_name,
324                        CS_FILE_MODE_WRITE,
325                        access[a_id],
326                        MPI_INFO_NULL,
327                        MPI_COMM_WORLD,
328                        MPI_COMM_WORLD);
329 
330 #else
331 
332       f = cs_file_open(output_file_name,
333                        CS_FILE_MODE_WRITE,
334                        access[a_id]);
335 
336 #endif /* (HAVE_MPI) */
337 
338       cs_file_set_big_endian(f);
339 
340       cs_file_dump(f);
341 
342       sprintf(buf, "fvm test file");
343       for (i = strlen(buf); i < 80; i++)
344         buf[i] = '\0';
345 
346       retval = cs_file_write_global(f, buf, 1, 80);
347 
348       if (rank == 0)
349         bft_printf("rank %d, wrote %d global values.\n", rank, (int)retval);
350 
351       for (i = block_start_2; i < block_end_2; i++) {
352         ibuf[(i-block_start_2)*2] = i*2 - 1;
353         ibuf[(i-block_start_2)*2 + 1] = i*2;
354       }
355       for (i = block_start; i < block_end; i++)
356         dbuf[i-block_start] = i;
357 
358       retval = cs_file_write_block(f, ibuf, sizeof(int), 2,
359                                    block_start_2, block_end_2);
360 
361       bft_printf("rank %d, wrote %d block values.\n", rank, (int)retval);
362 
363       retval = cs_file_write_block_buffer(f, dbuf, sizeof(double), 1,
364                                           block_start, block_end);
365 
366       bft_printf("rank %d, wrote %d block (buffer) values.\n",
367                  rank, (int)retval);
368 
369       sprintf(buf, "fvm test file end");
370       for (i = strlen(buf); i < 80; i++)
371         buf[i] = '\0';
372 
373       retval = cs_file_write_global(f, buf, 1, 80);
374 
375       if (rank == 0)
376         bft_printf("rank %d, wrote %d global values.\n", rank, (int)retval);
377 
378       f = cs_file_free(f);
379 
380       if (access[a_id] < CS_FILE_MPI_INDEPENDENT)
381         break;
382     }
383   }
384 
385   /* We are finished */
386 
387   bft_mem_end();
388 
389 #if defined(HAVE_MPI)
390 
391   MPI_Finalize();
392 
393 #endif
394 
395   exit (EXIT_SUCCESS);
396 }
397