xref: /phasta/phSolver/common/test/phIOreadIlwork.cc (revision 4d60bba28c1e1f3ca80b42756ae9dcbcd5c4bc48)
1 #include <mpi.h>
2 #include <stdio.h>
3 #include <stdlib.h>
4 #include <unistd.h>
5 #include <set>
6 #include "phIO.h"
7 #include "posixio.h"
8 #include "phio_posix.h"
9 
10 int main(int argc, char* argv[]) {
11   MPI_Init(&argc,&argv);
12   if( argc != 2 ) {
13     fprintf(stderr, "Usage: %s <geombc posix file>\n",argv[0]);
14     MPI_Finalize();
15     return 1;
16   }
17   const char* filename = argv[1];
18   const char* phrase = "ilwork";
19   const char* type = "integer";
20   const char* iotype = "binary";
21   int* ilwork = NULL;
22   int len = 0;
23   int one = 2;
24 
25   phio_fp file;
26   posixio_setup(&(file), 'r');
27   posix_openfile_single(filename, file);
28   phio_readheader(file, phrase, &len, &one, type, iotype);
29   fprintf(stderr, "len %d\n", len);
30   ilwork = (int*) malloc(len*sizeof(int));
31   phio_readdatablock(file, phrase, ilwork, &len, type, iotype);
32   phio_closefile(file);
33 
34   // Initialization
35   int itkbeg = 0;
36   int m = 0;
37   int idl = 0;
38   std::set<int> neighbors;
39 
40   // Compute summary info first such as number of communication tasks, neighboring parts, onwer parts, etc
41   int numtask = ilwork[0];
42   int numowner = 0;
43   for(int itask=0;itask<numtask;itask++) {
44     int iacc   = ilwork [itkbeg + 2];
45     int iother = ilwork [itkbeg + 3];
46     int numseg = ilwork [itkbeg + 4];
47     if(iacc == 1) numowner++;
48     neighbors.insert(iother);
49     itkbeg = itkbeg + 4 + 2*numseg;
50   }
51   printf("Number of neighboring parts: %d\n",neighbors.size());
52   printf("Number of communication tasks: %d\n",numtask);
53   printf("Number of owner communications: %d\n",numowner);
54   printf("Number of non-owner communications: %d\n",numtask-numowner);
55 
56   // Print now communication info
57   printf("\n");
58   printf("Communication info for this part:\n");
59   itkbeg = 0;
60   for(int itask=0;itask<numtask;itask++) {
61     int itag   = ilwork [itkbeg + 1];
62     int iacc   = ilwork [itkbeg + 2];
63     int iother = ilwork [itkbeg + 3];
64     int numseg = ilwork [itkbeg + 4];
65     int isgbeg = ilwork [itkbeg + 5];
66 
67     // Toal number of nodes involved in this communication (lfront), including all segments
68     int lfront = 0;
69     for(int is=1;is<=numseg;is++) {
70       int lenseg = ilwork [itkbeg + 4 + 2*is];
71       lfront = lfront + lenseg;
72     }
73     printf("Communication %d:\ttag %d\townership %d\trank %d\tnumseg %d\tnumvtx %d\n",itask,itag,iacc,iother,numseg,lfront);
74     itkbeg = itkbeg + 4 + 2*numseg;
75   }
76 
77   // Print now the raw ilwork array
78   printf("\n");
79   printf("ilwork array:\n");
80   for(int i=0;i<len;i++) {
81     printf("%d\n",ilwork[i]);
82   }
83 
84   free(ilwork);
85 
86   MPI_Finalize();
87   return 0;
88 }
89