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