xref: /petsc/src/sys/tests/ex17.c (revision 609caa7c8c030312b00807b4f015fd827bb80932)
1e9c46472SJames Wright static char help[] = "Demonstrates PetscFOpens() and PetscSynchronizedFGets().\n\n";
2e9c46472SJames Wright 
3e9c46472SJames Wright #include <petscsys.h>
main(int argc,char ** argv)4e9c46472SJames Wright int main(int argc, char **argv)
5e9c46472SJames Wright {
6e9c46472SJames Wright   const char  line1[]    = "hello 1\n";
7e9c46472SJames Wright   const char  line2[]    = "hello 2\n";
8e9c46472SJames Wright   const char  filename[] = "testfile";
9e9c46472SJames Wright   PetscMPIInt rank;
10e9c46472SJames Wright   FILE       *fp;
11e9c46472SJames Wright 
12e9c46472SJames Wright   PetscFunctionBeginUser;
13c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
14e9c46472SJames Wright   MPI_Comm comm = PETSC_COMM_WORLD;
15e9c46472SJames Wright   PetscCallMPI(MPI_Comm_rank(comm, &rank));
16e9c46472SJames Wright 
17e9c46472SJames Wright   // -- Create the file
18e9c46472SJames Wright   PetscCall(PetscFOpen(comm, filename, "w", &fp));
19e9c46472SJames Wright   PetscCall(PetscFPrintf(comm, fp, line1));
20e9c46472SJames Wright   PetscCall(PetscFPrintf(comm, fp, line2));
21e9c46472SJames Wright   PetscCall(PetscSynchronizedFPrintf(comm, fp, "rank: %d\n", rank)); // Print rankid in order
22e9c46472SJames Wright   PetscCall(PetscSynchronizedFlush(comm, fp));
23e9c46472SJames Wright   PetscCall(PetscFClose(comm, fp));
24e9c46472SJames Wright 
25e9c46472SJames Wright   { // -- Read the file
26e9c46472SJames Wright     char        line[512] = {0};
27e9c46472SJames Wright     PetscBool   line_check;
28e9c46472SJames Wright     PetscMPIInt size;
29e9c46472SJames Wright     PetscInt    line_rank;
30e9c46472SJames Wright 
31e9c46472SJames Wright     PetscCall(PetscFOpen(comm, filename, "r", &fp));
32e9c46472SJames Wright     PetscCall(PetscSynchronizedFGets(comm, fp, sizeof(line), line));
33e9c46472SJames Wright     PetscCall(PetscStrncmp(line, line1, sizeof(line1), &line_check));
34e9c46472SJames Wright     PetscCheck(line_check, PETSC_COMM_SELF, PETSC_ERR_FILE_READ, "Line 1 not read correctly. Got '%s', expected '%s'", line, line1);
35e9c46472SJames Wright     PetscCall(PetscSynchronizedFGets(comm, fp, sizeof(line), line));
36e9c46472SJames Wright     PetscCall(PetscStrncmp(line, line2, sizeof(line2), &line_check));
37e9c46472SJames Wright     PetscCheck(line_check, PETSC_COMM_SELF, PETSC_ERR_FILE_READ, "Line 2 not read correctly. Got '%s', expected '%s'", line, line2);
38e9c46472SJames Wright     PetscCallMPI(MPI_Comm_size(comm, &size));
39e9c46472SJames Wright     for (PetscInt i = 0; i < size; i++) {
40e9c46472SJames Wright       PetscCall(PetscSynchronizedFGets(comm, fp, sizeof(line), line));
41e9c46472SJames Wright       sscanf(line, "rank: %" PetscInt_FMT, &line_rank);
42e9c46472SJames Wright       PetscCheck(i == line_rank, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Did not find correct rank line in file. Expected %" PetscInt_FMT ", found %" PetscInt_FMT, i, line_rank);
43e9c46472SJames Wright     }
44e9c46472SJames Wright 
45e9c46472SJames Wright     PetscCall(PetscFClose(comm, fp));
46e9c46472SJames Wright   }
47e9c46472SJames Wright 
48e9c46472SJames Wright   PetscCall(PetscFinalize());
49e9c46472SJames Wright   return 0;
50e9c46472SJames Wright }
51e9c46472SJames Wright 
52e9c46472SJames Wright /*TEST
53e9c46472SJames Wright 
54e9c46472SJames Wright    test:
55e9c46472SJames Wright     nsize: 3
56*3886731fSPierre Jolivet     output_file: output/empty.out
57e9c46472SJames Wright 
58e9c46472SJames Wright TEST*/
59