1 static char help[] = "Test MatCreateSubMatrices\n\n"; 2 3 #include <petscis.h> 4 #include <petscmat.h> 5 6 int main(int argc,char **args) 7 { 8 PetscErrorCode ierr; 9 Mat A,*submats,*submats2; 10 IS *irow,*icol; 11 PetscInt i,n; 12 PetscMPIInt rank; 13 PetscViewer matfd,rowfd,colfd; 14 PetscBool same; 15 char matfile[PETSC_MAX_PATH_LEN],rowfile[PETSC_MAX_PATH_LEN],colfile[PETSC_MAX_PATH_LEN]; 16 char rankstr[16]={0}; 17 18 ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 19 CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 20 21 CHKERRQ(PetscOptionsGetString(NULL,NULL,"-A",matfile,sizeof(matfile),NULL)); 22 CHKERRQ(PetscOptionsGetString(NULL,NULL,"-row",rowfile,sizeof(rowfile),NULL)); 23 CHKERRQ(PetscOptionsGetString(NULL,NULL,"-col",colfile,sizeof(colfile),NULL)); 24 25 /* Each rank has its own files for row/col ISes */ 26 CHKERRQ(PetscSNPrintf(rankstr,16,"-%d",rank)); 27 CHKERRQ(PetscStrlcat(rowfile,rankstr,PETSC_MAX_PATH_LEN)); 28 CHKERRQ(PetscStrlcat(colfile,rankstr,PETSC_MAX_PATH_LEN)); 29 30 CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,matfile,FILE_MODE_READ,&matfd)); 31 CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_SELF,rowfile,FILE_MODE_READ,&rowfd)); 32 CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_SELF,colfile,FILE_MODE_READ,&colfd)); 33 34 CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 35 CHKERRQ(MatSetFromOptions(A)); 36 CHKERRQ(MatLoad(A,matfd)); 37 38 /* We stored the number of ISes at the beginning of rowfd */ 39 CHKERRQ(PetscViewerBinaryRead(rowfd,&n,1,NULL,PETSC_INT)); 40 CHKERRQ(PetscMalloc2(n,&irow,n,&icol)); 41 for (i=0; i<n; i++) { 42 CHKERRQ(ISCreate(PETSC_COMM_SELF,&irow[i])); 43 CHKERRQ(ISCreate(PETSC_COMM_SELF,&icol[i])); 44 CHKERRQ(ISLoad(irow[i],rowfd)); 45 CHKERRQ(ISLoad(icol[i],colfd)); 46 } 47 48 CHKERRQ(PetscViewerDestroy(&matfd)); 49 CHKERRQ(PetscViewerDestroy(&rowfd)); 50 CHKERRQ(PetscViewerDestroy(&colfd)); 51 52 /* Create submats for the first time */ 53 CHKERRQ(MatCreateSubMatrices(A,n,irow,icol,MAT_INITIAL_MATRIX,&submats)); 54 55 /* Dup submats to submats2 for later comparison */ 56 CHKERRQ(PetscMalloc1(n,&submats2)); 57 for (i=0; i<n; i++) { 58 CHKERRQ(MatDuplicate(submats[i],MAT_COPY_VALUES,&submats2[i])); 59 } 60 61 /* Create submats again */ 62 CHKERRQ(MatCreateSubMatrices(A,n,irow,icol,MAT_REUSE_MATRIX,&submats)); 63 64 /* Compare submats and submats2 */ 65 for (i=0; i<n; i++) { 66 CHKERRQ(MatEqual(submats[i],submats2[i],&same)); 67 PetscCheck(same,PETSC_COMM_SELF,PETSC_ERR_PLIB,"submatrix %" PetscInt_FMT " is not same",i); 68 } 69 70 CHKERRQ(MatDestroy(&A)); 71 for (i=0; i<n; i++) { 72 CHKERRQ(ISDestroy(&irow[i])); 73 CHKERRQ(ISDestroy(&icol[i])); 74 } 75 CHKERRQ(MatDestroySubMatrices(n,&submats)); 76 CHKERRQ(MatDestroyMatrices(n,&submats2)); 77 CHKERRQ(PetscFree2(irow,icol)); 78 ierr = PetscFinalize(); 79 return ierr; 80 } 81 82 /*TEST 83 84 test: 85 suffix: 1 86 nsize: 2 87 requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) 88 args: -mat_type {{aij baij}} -A ${DATAFILESPATH}/matrices/CreateSubMatrices/A -row ${DATAFILESPATH}/matrices/CreateSubMatrices/row -col ${DATAFILESPATH}/matrices/CreateSubMatrices/col 89 90 TEST*/ 91