1544b287aSSatish Balay static char help[] = "Test MatCreateSubMatrices\n\n"; 2544b287aSSatish Balay 3544b287aSSatish Balay #include <petscis.h> 4544b287aSSatish Balay #include <petscmat.h> 5544b287aSSatish Balay 6544b287aSSatish Balay int main(int argc,char **args) 7544b287aSSatish Balay { 8544b287aSSatish Balay Mat A,*submats,*submats2; 9544b287aSSatish Balay IS *irow,*icol; 10544b287aSSatish Balay PetscInt i,n; 11544b287aSSatish Balay PetscMPIInt rank; 12544b287aSSatish Balay PetscViewer matfd,rowfd,colfd; 13544b287aSSatish Balay PetscBool same; 14544b287aSSatish Balay char matfile[PETSC_MAX_PATH_LEN],rowfile[PETSC_MAX_PATH_LEN],colfile[PETSC_MAX_PATH_LEN]; 15544b287aSSatish Balay char rankstr[16]={0}; 16544b287aSSatish Balay 17*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 18*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 19544b287aSSatish Balay 20*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL,NULL,"-A",matfile,sizeof(matfile),NULL)); 21*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL,NULL,"-row",rowfile,sizeof(rowfile),NULL)); 22*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL,NULL,"-col",colfile,sizeof(colfile),NULL)); 23544b287aSSatish Balay 24544b287aSSatish Balay /* Each rank has its own files for row/col ISes */ 25*9566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(rankstr,16,"-%d",rank)); 26*9566063dSJacob Faibussowitsch PetscCall(PetscStrlcat(rowfile,rankstr,PETSC_MAX_PATH_LEN)); 27*9566063dSJacob Faibussowitsch PetscCall(PetscStrlcat(colfile,rankstr,PETSC_MAX_PATH_LEN)); 28544b287aSSatish Balay 29*9566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,matfile,FILE_MODE_READ,&matfd)); 30*9566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_SELF,rowfile,FILE_MODE_READ,&rowfd)); 31*9566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_SELF,colfile,FILE_MODE_READ,&colfd)); 32544b287aSSatish Balay 33*9566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); 34*9566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 35*9566063dSJacob Faibussowitsch PetscCall(MatLoad(A,matfd)); 36544b287aSSatish Balay 37544b287aSSatish Balay /* We stored the number of ISes at the beginning of rowfd */ 38*9566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(rowfd,&n,1,NULL,PETSC_INT)); 39*9566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(n,&irow,n,&icol)); 40544b287aSSatish Balay for (i=0; i<n; i++) { 41*9566063dSJacob Faibussowitsch PetscCall(ISCreate(PETSC_COMM_SELF,&irow[i])); 42*9566063dSJacob Faibussowitsch PetscCall(ISCreate(PETSC_COMM_SELF,&icol[i])); 43*9566063dSJacob Faibussowitsch PetscCall(ISLoad(irow[i],rowfd)); 44*9566063dSJacob Faibussowitsch PetscCall(ISLoad(icol[i],colfd)); 45544b287aSSatish Balay } 46544b287aSSatish Balay 47*9566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&matfd)); 48*9566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&rowfd)); 49*9566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&colfd)); 50544b287aSSatish Balay 51544b287aSSatish Balay /* Create submats for the first time */ 52*9566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(A,n,irow,icol,MAT_INITIAL_MATRIX,&submats)); 53544b287aSSatish Balay 54544b287aSSatish Balay /* Dup submats to submats2 for later comparison */ 55*9566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n,&submats2)); 56544b287aSSatish Balay for (i=0; i<n; i++) { 57*9566063dSJacob Faibussowitsch PetscCall(MatDuplicate(submats[i],MAT_COPY_VALUES,&submats2[i])); 58544b287aSSatish Balay } 59544b287aSSatish Balay 60544b287aSSatish Balay /* Create submats again */ 61*9566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(A,n,irow,icol,MAT_REUSE_MATRIX,&submats)); 62544b287aSSatish Balay 63544b287aSSatish Balay /* Compare submats and submats2 */ 64544b287aSSatish Balay for (i=0; i<n; i++) { 65*9566063dSJacob Faibussowitsch PetscCall(MatEqual(submats[i],submats2[i],&same)); 6628b400f6SJacob Faibussowitsch PetscCheck(same,PETSC_COMM_SELF,PETSC_ERR_PLIB,"submatrix %" PetscInt_FMT " is not same",i); 67544b287aSSatish Balay } 68544b287aSSatish Balay 69*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 70544b287aSSatish Balay for (i=0; i<n; i++) { 71*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(&irow[i])); 72*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(&icol[i])); 73544b287aSSatish Balay } 74*9566063dSJacob Faibussowitsch PetscCall(MatDestroySubMatrices(n,&submats)); 75*9566063dSJacob Faibussowitsch PetscCall(MatDestroyMatrices(n,&submats2)); 76*9566063dSJacob Faibussowitsch PetscCall(PetscFree2(irow,icol)); 77*9566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 78b122ec5aSJacob Faibussowitsch return 0; 79544b287aSSatish Balay } 80544b287aSSatish Balay 81544b287aSSatish Balay /*TEST 82544b287aSSatish Balay 83544b287aSSatish Balay test: 84544b287aSSatish Balay suffix: 1 85544b287aSSatish Balay nsize: 2 86dfd57a17SPierre Jolivet requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) 87544b287aSSatish Balay args: -mat_type {{aij baij}} -A ${DATAFILESPATH}/matrices/CreateSubMatrices/A -row ${DATAFILESPATH}/matrices/CreateSubMatrices/row -col ${DATAFILESPATH}/matrices/CreateSubMatrices/col 88544b287aSSatish Balay 89544b287aSSatish Balay TEST*/ 90