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