xref: /phasta/phSolver/common/test/phIOreadIlwork.cc (revision a58fb009befd71b3764bdd319273b61df5485cf0)
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 
13560e081fSCameron Smith int main(int argc, char* argv[]) {
14560e081fSCameron Smith   MPI_Init(&argc,&argv);
15302b84e8SCameron Smith   int commrank,commsize;
16302b84e8SCameron Smith   MPI_Comm_rank(MPI_COMM_WORLD,&commrank);
17302b84e8SCameron Smith   MPI_Comm_size(MPI_COMM_WORLD,&commsize);
18eeef5fd6SCameron Smith   if( argc != 5 ) {
19302b84e8SCameron Smith     if( !commrank )
20eeef5fd6SCameron Smith       fprintf(stderr, "Usage: %s <geombc posix file> <verbosity=0|1|2> <rankoffset> <outfile>\n",argv[0]);
21eeef5fd6SCameron Smith       fprintf(stderr, "verbosity=0 will only write to the specified \'outfile\'\n",argv[0]);
22eeef5fd6SCameron Smith       fprintf(stderr, "verbosity>0 will write to the specified \'outfile\' and to stdout\n",argv[0]);
23560e081fSCameron Smith     MPI_Finalize();
24560e081fSCameron Smith     return 1;
25560e081fSCameron Smith   }
26560e081fSCameron Smith   const char* filename = argv[1];
27302b84e8SCameron Smith   const int verbosity = atoi(argv[2]);
28302b84e8SCameron Smith   const int rankoffset = atoi(argv[3]);
29eeef5fd6SCameron Smith   char* outfilename = argv[4];
30560e081fSCameron Smith   const char* phrase = "ilwork";
31560e081fSCameron Smith   const char* type = "integer";
32560e081fSCameron Smith   const char* iotype = "binary";
33560e081fSCameron Smith   int* ilwork = NULL;
34560e081fSCameron Smith   int len = 0;
35560e081fSCameron Smith   int one = 2;
36560e081fSCameron Smith 
37eb8873f5SCameron Smith   const int phastarank = commrank+1+rankoffset;
38eb8873f5SCameron Smith   const int group = (commrank+rankoffset)/2048;
39302b84e8SCameron Smith 
40302b84e8SCameron Smith   char geomfilename[4096];
41eb8873f5SCameron Smith   sprintf(geomfilename, "%s/%d/geombc.dat.%d",filename,group,phastarank);
42302b84e8SCameron Smith 
43560e081fSCameron Smith   phio_fp file;
44560e081fSCameron Smith   posixio_setup(&(file), 'r');
45302b84e8SCameron Smith   posix_openfile_single(geomfilename, file);
46560e081fSCameron Smith   phio_readheader(file, phrase, &len, &one, type, iotype);
47560e081fSCameron Smith   ilwork = (int*) malloc(len*sizeof(int));
48560e081fSCameron Smith   phio_readdatablock(file, phrase, ilwork, &len, type, iotype);
49560e081fSCameron Smith   phio_closefile(file);
504894d94eSMichel Rasquin 
514894d94eSMichel Rasquin   // Initialization
524894d94eSMichel Rasquin   int itkbeg = 0;
534894d94eSMichel Rasquin   int m = 0;
544894d94eSMichel Rasquin   int idl = 0;
554894d94eSMichel Rasquin   std::set<int> neighbors;
564894d94eSMichel Rasquin 
574894d94eSMichel Rasquin   // Compute summary info first such as number of communication tasks, neighboring parts, onwer parts, etc
584894d94eSMichel Rasquin   int numtask = ilwork[0];
594894d94eSMichel Rasquin   int numowner = 0;
604894d94eSMichel Rasquin   for(int itask=0;itask<numtask;itask++) {
614894d94eSMichel Rasquin     int iacc   = ilwork [itkbeg + 2];
624894d94eSMichel Rasquin     int iother = ilwork [itkbeg + 3];
634894d94eSMichel Rasquin     int numseg = ilwork [itkbeg + 4];
644894d94eSMichel Rasquin     if(iacc == 1) numowner++;
654894d94eSMichel Rasquin     neighbors.insert(iother);
664894d94eSMichel Rasquin     itkbeg = itkbeg + 4 + 2*numseg;
674894d94eSMichel Rasquin   }
68c0903dacSCameron Smith   assert(neighbors.size() != 0);
69eeef5fd6SCameron Smith   MPI_Status status;
70eeef5fd6SCameron Smith   MPI_File outfile;
71eeef5fd6SCameron Smith   MPI_File_open(MPI_COMM_WORLD,outfilename,
72eeef5fd6SCameron Smith       MPI_MODE_CREATE|MPI_MODE_WRONLY,MPI_INFO_NULL,&outfile);
73c0903dacSCameron Smith   std::string header("rank,peers,tasks,owned,notowned\n");
74eeef5fd6SCameron Smith   if( !commrank ) //write header
75*a58fb009SCameron Smith     MPI_File_write_at(outfile,0,(void*)header.c_str(),header.size(),MPI_BYTE,&status);
76eeef5fd6SCameron Smith   std::stringstream ss;
77eeef5fd6SCameron Smith   ss << phastarank << ","
78eeef5fd6SCameron Smith      << neighbors.size() << ","
79eeef5fd6SCameron Smith      << numtask << ","
80eeef5fd6SCameron Smith      << numowner << ","
81eeef5fd6SCameron Smith      << numtask-numowner << "\n";
82eeef5fd6SCameron Smith   std::string s = ss.str();
83eeef5fd6SCameron Smith   int size = s.size();
84eeef5fd6SCameron Smith   int offset = 0;
85eeef5fd6SCameron Smith   MPI_Exscan(&size,&offset,1,MPI_INT,MPI_SUM,MPI_COMM_WORLD);
86eeef5fd6SCameron Smith   offset += header.size();
87*a58fb009SCameron Smith   int ret = MPI_File_write_at(outfile,offset,(void*)s.c_str(),s.size(),MPI_BYTE,&status);
88c0903dacSCameron Smith   assert(ret == MPI_SUCCESS);
89302b84e8SCameron Smith   if( verbosity > 0 ) {
904894d94eSMichel Rasquin     // Print now communication info
914894d94eSMichel Rasquin     printf("\n");
924894d94eSMichel Rasquin     printf("Communication info for this part:\n");
934894d94eSMichel Rasquin     itkbeg = 0;
944894d94eSMichel Rasquin     for(int itask=0;itask<numtask;itask++) {
954894d94eSMichel Rasquin       int itag   = ilwork [itkbeg + 1];
964894d94eSMichel Rasquin       int iacc   = ilwork [itkbeg + 2];
974894d94eSMichel Rasquin       int iother = ilwork [itkbeg + 3];
984894d94eSMichel Rasquin       int numseg = ilwork [itkbeg + 4];
994894d94eSMichel Rasquin       int isgbeg = ilwork [itkbeg + 5];
1004894d94eSMichel Rasquin 
1014894d94eSMichel Rasquin       // Toal number of nodes involved in this communication (lfront), including all segments
1024894d94eSMichel Rasquin       int lfront = 0;
1034894d94eSMichel Rasquin       for(int is=1;is<=numseg;is++) {
1044894d94eSMichel Rasquin         int lenseg = ilwork [itkbeg + 4 + 2*is];
1054894d94eSMichel Rasquin         lfront = lfront + lenseg;
1064894d94eSMichel Rasquin       }
1074894d94eSMichel Rasquin       printf("Communication %d:\ttag %d\townership %d\trank %d\tnumseg %d\tnumvtx %d\n",itask,itag,iacc,iother,numseg,lfront);
1084894d94eSMichel Rasquin       itkbeg = itkbeg + 4 + 2*numseg;
1094894d94eSMichel Rasquin     }
110302b84e8SCameron Smith     printf("\n");
111302b84e8SCameron Smith   }
1124894d94eSMichel Rasquin 
1134894d94eSMichel Rasquin   // Print now the raw ilwork array
114302b84e8SCameron Smith   if( verbosity > 1 ) {
1154894d94eSMichel Rasquin     printf("ilwork array:\n");
1164894d94eSMichel Rasquin     for(int i=0;i<len;i++) {
1174894d94eSMichel Rasquin       printf("%d\n",ilwork[i]);
1184894d94eSMichel Rasquin     }
11925b215f0SCameron Smith   }
1204894d94eSMichel Rasquin 
121560e081fSCameron Smith   free(ilwork);
122eeef5fd6SCameron Smith   MPI_File_close(&outfile);
123560e081fSCameron Smith   MPI_Finalize();
124560e081fSCameron Smith   return 0;
125560e081fSCameron Smith }
126