xref: /phasta/phSolver/common/test/phIOreadIlwork.cc (revision 021b6b1dd9f35c31a2282898b06fbc7e3008e3dc)
1560e081fSCameron Smith #include <mpi.h>
2560e081fSCameron Smith #include <stdio.h>
3560e081fSCameron Smith #include <stdlib.h>
4560e081fSCameron Smith #include <unistd.h>
54894d94eSMichel Rasquin #include <set>
6eeef5fd6SCameron Smith #include <string>
7eeef5fd6SCameron Smith #include <sstream>
8c0903dacSCameron Smith #include <cassert>
9560e081fSCameron Smith #include "phIO.h"
10560e081fSCameron Smith #include "posixio.h"
11560e081fSCameron Smith #include "phio_posix.h"
12560e081fSCameron Smith 
openfile(const char * geomfilename,phio_fp & file)13ff414521SCameron Smith int openfile(const char* geomfilename, phio_fp& file) {
14ff414521SCameron Smith   int commrank;
15ff414521SCameron Smith   MPI_Comm_rank(MPI_COMM_WORLD,&commrank);
16ff414521SCameron Smith   int err = posix_openfile_single(geomfilename, file);
17ff414521SCameron Smith   int globalerr = 0;
18ff414521SCameron Smith   MPI_Allreduce(&err, &globalerr, 1, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
19ff414521SCameron Smith   if(err > 0 && !commrank)
20ff414521SCameron Smith     fprintf(stderr, "failed to open file %s\n", geomfilename);
21ff414521SCameron Smith   return globalerr;
22ff414521SCameron Smith }
23ff414521SCameron Smith 
main(int argc,char * argv[])24560e081fSCameron Smith int main(int argc, char* argv[]) {
25560e081fSCameron Smith   MPI_Init(&argc,&argv);
26302b84e8SCameron Smith   int commrank,commsize;
27302b84e8SCameron Smith   MPI_Comm_rank(MPI_COMM_WORLD,&commrank);
28302b84e8SCameron Smith   MPI_Comm_size(MPI_COMM_WORLD,&commsize);
29eeef5fd6SCameron Smith   if( argc != 5 ) {
30302b84e8SCameron Smith     if( !commrank )
31eeef5fd6SCameron Smith       fprintf(stderr, "Usage: %s <geombc posix file> <verbosity=0|1|2> <rankoffset> <outfile>\n",argv[0]);
32*021b6b1dSCameron Smith       fprintf(stderr, "verbosity=0 will only write to the specified \'outfile\'\n");
33*021b6b1dSCameron Smith       fprintf(stderr, "verbosity>0 will write to the specified \'outfile\' and to stdout\n");
34560e081fSCameron Smith     MPI_Finalize();
35560e081fSCameron Smith     return 1;
36560e081fSCameron Smith   }
37560e081fSCameron Smith   const char* filename = argv[1];
38302b84e8SCameron Smith   const int verbosity = atoi(argv[2]);
39302b84e8SCameron Smith   const int rankoffset = atoi(argv[3]);
40eeef5fd6SCameron Smith   char* outfilename = argv[4];
41560e081fSCameron Smith   const char* phrase = "ilwork";
42560e081fSCameron Smith   const char* type = "integer";
43560e081fSCameron Smith   const char* iotype = "binary";
44560e081fSCameron Smith   int* ilwork = NULL;
45560e081fSCameron Smith   int len = 0;
46560e081fSCameron Smith   int one = 2;
47560e081fSCameron Smith 
48eb8873f5SCameron Smith   const int phastarank = commrank+1+rankoffset;
49eb8873f5SCameron Smith   const int group = (commrank+rankoffset)/2048;
50302b84e8SCameron Smith 
51302b84e8SCameron Smith   char geomfilename[4096];
52eb8873f5SCameron Smith   sprintf(geomfilename, "%s/%d/geombc.dat.%d",filename,group,phastarank);
53302b84e8SCameron Smith 
54560e081fSCameron Smith   phio_fp file;
55560e081fSCameron Smith   posixio_setup(&(file), 'r');
56ff414521SCameron Smith   int err = openfile(geomfilename, file);
57ff414521SCameron Smith   if(err > 0) {
58ff414521SCameron Smith     if(!commrank) fprintf(stderr, "trying again without the sub-directory...\n");
59ff414521SCameron Smith     sprintf(geomfilename, "%s/geombc.dat.%d",filename,phastarank);
60ff414521SCameron Smith     err = openfile(geomfilename, file);
61ff414521SCameron Smith     assert(!err);
62ff414521SCameron Smith     if(!commrank)
63ff414521SCameron Smith       fprintf(stderr, "geombc files opened successfully\n");
64ff414521SCameron Smith   }
65560e081fSCameron Smith   phio_readheader(file, phrase, &len, &one, type, iotype);
66560e081fSCameron Smith   ilwork = (int*) malloc(len*sizeof(int));
67560e081fSCameron Smith   phio_readdatablock(file, phrase, ilwork, &len, type, iotype);
68560e081fSCameron Smith   phio_closefile(file);
694894d94eSMichel Rasquin 
704894d94eSMichel Rasquin   // Initialization
714894d94eSMichel Rasquin   int itkbeg = 0;
724894d94eSMichel Rasquin   std::set<int> neighbors;
734894d94eSMichel Rasquin 
744894d94eSMichel Rasquin   // Compute summary info first such as number of communication tasks, neighboring parts, onwer parts, etc
754894d94eSMichel Rasquin   int numtask = ilwork[0];
764894d94eSMichel Rasquin   int numowner = 0;
774894d94eSMichel Rasquin   for(int itask=0;itask<numtask;itask++) {
784894d94eSMichel Rasquin     int iacc   = ilwork [itkbeg + 2];
794894d94eSMichel Rasquin     int iother = ilwork [itkbeg + 3];
804894d94eSMichel Rasquin     int numseg = ilwork [itkbeg + 4];
814894d94eSMichel Rasquin     if(iacc == 1) numowner++;
824894d94eSMichel Rasquin     neighbors.insert(iother);
834894d94eSMichel Rasquin     itkbeg = itkbeg + 4 + 2*numseg;
844894d94eSMichel Rasquin   }
85c0903dacSCameron Smith   assert(neighbors.size() != 0);
86eeef5fd6SCameron Smith   MPI_Status status;
87eeef5fd6SCameron Smith   MPI_File outfile;
88eeef5fd6SCameron Smith   MPI_File_open(MPI_COMM_WORLD,outfilename,
89eeef5fd6SCameron Smith       MPI_MODE_CREATE|MPI_MODE_WRONLY,MPI_INFO_NULL,&outfile);
90c0903dacSCameron Smith   std::string header("rank,peers,tasks,owned,notowned\n");
91eeef5fd6SCameron Smith   if( !commrank ) //write header
92a58fb009SCameron Smith     MPI_File_write_at(outfile,0,(void*)header.c_str(),header.size(),MPI_BYTE,&status);
93eeef5fd6SCameron Smith   std::stringstream ss;
94eeef5fd6SCameron Smith   ss << phastarank << ","
95eeef5fd6SCameron Smith      << neighbors.size() << ","
96eeef5fd6SCameron Smith      << numtask << ","
97eeef5fd6SCameron Smith      << numowner << ","
98eeef5fd6SCameron Smith      << numtask-numowner << "\n";
99eeef5fd6SCameron Smith   std::string s = ss.str();
100eeef5fd6SCameron Smith   int size = s.size();
101eeef5fd6SCameron Smith   int offset = 0;
102eeef5fd6SCameron Smith   MPI_Exscan(&size,&offset,1,MPI_INT,MPI_SUM,MPI_COMM_WORLD);
103eeef5fd6SCameron Smith   offset += header.size();
104a58fb009SCameron Smith   int ret = MPI_File_write_at(outfile,offset,(void*)s.c_str(),s.size(),MPI_BYTE,&status);
105c0903dacSCameron Smith   assert(ret == MPI_SUCCESS);
106302b84e8SCameron Smith   if( verbosity > 0 ) {
1074894d94eSMichel Rasquin     // Print now communication info
1084894d94eSMichel Rasquin     printf("\n");
1094894d94eSMichel Rasquin     printf("Communication info for this part:\n");
1104894d94eSMichel Rasquin     itkbeg = 0;
1114894d94eSMichel Rasquin     for(int itask=0;itask<numtask;itask++) {
1124894d94eSMichel Rasquin       int itag   = ilwork [itkbeg + 1];
1134894d94eSMichel Rasquin       int iacc   = ilwork [itkbeg + 2];
1144894d94eSMichel Rasquin       int iother = ilwork [itkbeg + 3];
1154894d94eSMichel Rasquin       int numseg = ilwork [itkbeg + 4];
1164894d94eSMichel Rasquin 
1174894d94eSMichel Rasquin       // Toal number of nodes involved in this communication (lfront), including all segments
1184894d94eSMichel Rasquin       int lfront = 0;
1194894d94eSMichel Rasquin       for(int is=1;is<=numseg;is++) {
1204894d94eSMichel Rasquin         int lenseg = ilwork [itkbeg + 4 + 2*is];
1214894d94eSMichel Rasquin         lfront = lfront + lenseg;
1224894d94eSMichel Rasquin       }
1234894d94eSMichel Rasquin       printf("Communication %d:\ttag %d\townership %d\trank %d\tnumseg %d\tnumvtx %d\n",itask,itag,iacc,iother,numseg,lfront);
1244894d94eSMichel Rasquin       itkbeg = itkbeg + 4 + 2*numseg;
1254894d94eSMichel Rasquin     }
126302b84e8SCameron Smith     printf("\n");
127302b84e8SCameron Smith   }
1284894d94eSMichel Rasquin 
1294894d94eSMichel Rasquin   // Print now the raw ilwork array
130302b84e8SCameron Smith   if( verbosity > 1 ) {
1314894d94eSMichel Rasquin     printf("ilwork array:\n");
1324894d94eSMichel Rasquin     for(int i=0;i<len;i++) {
1334894d94eSMichel Rasquin       printf("%d\n",ilwork[i]);
1344894d94eSMichel Rasquin     }
13525b215f0SCameron Smith   }
1364894d94eSMichel Rasquin 
137560e081fSCameron Smith   free(ilwork);
138eeef5fd6SCameron Smith   MPI_File_close(&outfile);
139560e081fSCameron Smith   MPI_Finalize();
140560e081fSCameron Smith   return 0;
141560e081fSCameron Smith }
142