xref: /phasta/phSolver/common/test/phIOreadIlwork.cc (revision 021b6b1dd9f35c31a2282898b06fbc7e3008e3dc)
1 #include <mpi.h>
2 #include <stdio.h>
3 #include <stdlib.h>
4 #include <unistd.h>
5 #include <set>
6 #include <string>
7 #include <sstream>
8 #include <cassert>
9 #include "phIO.h"
10 #include "posixio.h"
11 #include "phio_posix.h"
12 
openfile(const char * geomfilename,phio_fp & file)13 int openfile(const char* geomfilename, phio_fp& file) {
14   int commrank;
15   MPI_Comm_rank(MPI_COMM_WORLD,&commrank);
16   int err = posix_openfile_single(geomfilename, file);
17   int globalerr = 0;
18   MPI_Allreduce(&err, &globalerr, 1, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
19   if(err > 0 && !commrank)
20     fprintf(stderr, "failed to open file %s\n", geomfilename);
21   return globalerr;
22 }
23 
main(int argc,char * argv[])24 int main(int argc, char* argv[]) {
25   MPI_Init(&argc,&argv);
26   int commrank,commsize;
27   MPI_Comm_rank(MPI_COMM_WORLD,&commrank);
28   MPI_Comm_size(MPI_COMM_WORLD,&commsize);
29   if( argc != 5 ) {
30     if( !commrank )
31       fprintf(stderr, "Usage: %s <geombc posix file> <verbosity=0|1|2> <rankoffset> <outfile>\n",argv[0]);
32       fprintf(stderr, "verbosity=0 will only write to the specified \'outfile\'\n");
33       fprintf(stderr, "verbosity>0 will write to the specified \'outfile\' and to stdout\n");
34     MPI_Finalize();
35     return 1;
36   }
37   const char* filename = argv[1];
38   const int verbosity = atoi(argv[2]);
39   const int rankoffset = atoi(argv[3]);
40   char* outfilename = argv[4];
41   const char* phrase = "ilwork";
42   const char* type = "integer";
43   const char* iotype = "binary";
44   int* ilwork = NULL;
45   int len = 0;
46   int one = 2;
47 
48   const int phastarank = commrank+1+rankoffset;
49   const int group = (commrank+rankoffset)/2048;
50 
51   char geomfilename[4096];
52   sprintf(geomfilename, "%s/%d/geombc.dat.%d",filename,group,phastarank);
53 
54   phio_fp file;
55   posixio_setup(&(file), 'r');
56   int err = openfile(geomfilename, file);
57   if(err > 0) {
58     if(!commrank) fprintf(stderr, "trying again without the sub-directory...\n");
59     sprintf(geomfilename, "%s/geombc.dat.%d",filename,phastarank);
60     err = openfile(geomfilename, file);
61     assert(!err);
62     if(!commrank)
63       fprintf(stderr, "geombc files opened successfully\n");
64   }
65   phio_readheader(file, phrase, &len, &one, type, iotype);
66   ilwork = (int*) malloc(len*sizeof(int));
67   phio_readdatablock(file, phrase, ilwork, &len, type, iotype);
68   phio_closefile(file);
69 
70   // Initialization
71   int itkbeg = 0;
72   std::set<int> neighbors;
73 
74   // Compute summary info first such as number of communication tasks, neighboring parts, onwer parts, etc
75   int numtask = ilwork[0];
76   int numowner = 0;
77   for(int itask=0;itask<numtask;itask++) {
78     int iacc   = ilwork [itkbeg + 2];
79     int iother = ilwork [itkbeg + 3];
80     int numseg = ilwork [itkbeg + 4];
81     if(iacc == 1) numowner++;
82     neighbors.insert(iother);
83     itkbeg = itkbeg + 4 + 2*numseg;
84   }
85   assert(neighbors.size() != 0);
86   MPI_Status status;
87   MPI_File outfile;
88   MPI_File_open(MPI_COMM_WORLD,outfilename,
89       MPI_MODE_CREATE|MPI_MODE_WRONLY,MPI_INFO_NULL,&outfile);
90   std::string header("rank,peers,tasks,owned,notowned\n");
91   if( !commrank ) //write header
92     MPI_File_write_at(outfile,0,(void*)header.c_str(),header.size(),MPI_BYTE,&status);
93   std::stringstream ss;
94   ss << phastarank << ","
95      << neighbors.size() << ","
96      << numtask << ","
97      << numowner << ","
98      << numtask-numowner << "\n";
99   std::string s = ss.str();
100   int size = s.size();
101   int offset = 0;
102   MPI_Exscan(&size,&offset,1,MPI_INT,MPI_SUM,MPI_COMM_WORLD);
103   offset += header.size();
104   int ret = MPI_File_write_at(outfile,offset,(void*)s.c_str(),s.size(),MPI_BYTE,&status);
105   assert(ret == MPI_SUCCESS);
106   if( verbosity > 0 ) {
107     // Print now communication info
108     printf("\n");
109     printf("Communication info for this part:\n");
110     itkbeg = 0;
111     for(int itask=0;itask<numtask;itask++) {
112       int itag   = ilwork [itkbeg + 1];
113       int iacc   = ilwork [itkbeg + 2];
114       int iother = ilwork [itkbeg + 3];
115       int numseg = ilwork [itkbeg + 4];
116 
117       // Toal number of nodes involved in this communication (lfront), including all segments
118       int lfront = 0;
119       for(int is=1;is<=numseg;is++) {
120         int lenseg = ilwork [itkbeg + 4 + 2*is];
121         lfront = lfront + lenseg;
122       }
123       printf("Communication %d:\ttag %d\townership %d\trank %d\tnumseg %d\tnumvtx %d\n",itask,itag,iacc,iother,numseg,lfront);
124       itkbeg = itkbeg + 4 + 2*numseg;
125     }
126     printf("\n");
127   }
128 
129   // Print now the raw ilwork array
130   if( verbosity > 1 ) {
131     printf("ilwork array:\n");
132     for(int i=0;i<len;i++) {
133       printf("%d\n",ilwork[i]);
134     }
135   }
136 
137   free(ilwork);
138   MPI_File_close(&outfile);
139   MPI_Finalize();
140   return 0;
141 }
142