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 PetscErrorCode ierr; 9544b287aSSatish Balay Mat A,*submats,*submats2; 10544b287aSSatish Balay IS *irow,*icol; 11544b287aSSatish Balay PetscInt i,n; 12544b287aSSatish Balay PetscMPIInt rank; 13544b287aSSatish Balay PetscViewer matfd,rowfd,colfd; 14544b287aSSatish Balay PetscBool same; 15544b287aSSatish Balay char matfile[PETSC_MAX_PATH_LEN],rowfile[PETSC_MAX_PATH_LEN],colfile[PETSC_MAX_PATH_LEN]; 16544b287aSSatish Balay char rankstr[16]={0}; 17544b287aSSatish Balay 18544b287aSSatish Balay ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 19544b287aSSatish Balay ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr); 20544b287aSSatish Balay 21544b287aSSatish Balay ierr = PetscOptionsGetString(NULL,NULL,"-A",matfile,sizeof(matfile),NULL);CHKERRQ(ierr); 22544b287aSSatish Balay ierr = PetscOptionsGetString(NULL,NULL,"-row",rowfile,sizeof(rowfile),NULL);CHKERRQ(ierr); 23544b287aSSatish Balay ierr = PetscOptionsGetString(NULL,NULL,"-col",colfile,sizeof(colfile),NULL);CHKERRQ(ierr); 24544b287aSSatish Balay 25544b287aSSatish Balay /* Each rank has its own files for row/col ISes */ 26544b287aSSatish Balay ierr = PetscSNPrintf(rankstr,16,"-%d",rank);CHKERRQ(ierr); 27544b287aSSatish Balay ierr = PetscStrlcat(rowfile,rankstr,PETSC_MAX_PATH_LEN);CHKERRQ(ierr); 28544b287aSSatish Balay ierr = PetscStrlcat(colfile,rankstr,PETSC_MAX_PATH_LEN);CHKERRQ(ierr); 29544b287aSSatish Balay 30544b287aSSatish Balay ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,matfile,FILE_MODE_READ,&matfd);CHKERRQ(ierr); 31544b287aSSatish Balay ierr = PetscViewerBinaryOpen(PETSC_COMM_SELF,rowfile,FILE_MODE_READ,&rowfd);CHKERRQ(ierr); 32544b287aSSatish Balay ierr = PetscViewerBinaryOpen(PETSC_COMM_SELF,colfile,FILE_MODE_READ,&colfd);CHKERRQ(ierr); 33544b287aSSatish Balay 34544b287aSSatish Balay ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 35544b287aSSatish Balay ierr = MatSetFromOptions(A);CHKERRQ(ierr); 36544b287aSSatish Balay ierr = MatLoad(A,matfd);CHKERRQ(ierr); 37544b287aSSatish Balay 38544b287aSSatish Balay /* We stored the number of ISes at the beginning of rowfd */ 39544b287aSSatish Balay ierr = PetscViewerBinaryRead(rowfd,&n,1,NULL,PETSC_INT);CHKERRQ(ierr); 40544b287aSSatish Balay ierr = PetscMalloc2(n,&irow,n,&icol);CHKERRQ(ierr); 41544b287aSSatish Balay for (i=0; i<n; i++) { 42544b287aSSatish Balay ierr = ISCreate(PETSC_COMM_SELF,&irow[i]);CHKERRQ(ierr); 43544b287aSSatish Balay ierr = ISCreate(PETSC_COMM_SELF,&icol[i]);CHKERRQ(ierr); 44544b287aSSatish Balay ierr = ISLoad(irow[i],rowfd);CHKERRQ(ierr); 45544b287aSSatish Balay ierr = ISLoad(icol[i],colfd);CHKERRQ(ierr); 46544b287aSSatish Balay } 47544b287aSSatish Balay 48544b287aSSatish Balay ierr = PetscViewerDestroy(&matfd);CHKERRQ(ierr); 49544b287aSSatish Balay ierr = PetscViewerDestroy(&rowfd);CHKERRQ(ierr); 50544b287aSSatish Balay ierr = PetscViewerDestroy(&colfd);CHKERRQ(ierr); 51544b287aSSatish Balay 52544b287aSSatish Balay /* Create submats for the first time */ 53544b287aSSatish Balay ierr = MatCreateSubMatrices(A,n,irow,icol,MAT_INITIAL_MATRIX,&submats);CHKERRQ(ierr); 54544b287aSSatish Balay 55544b287aSSatish Balay /* Dup submats to submats2 for later comparison */ 56544b287aSSatish Balay ierr = PetscMalloc1(n,&submats2);CHKERRQ(ierr); 57544b287aSSatish Balay for (i=0; i<n; i++) { 58544b287aSSatish Balay ierr = MatDuplicate(submats[i],MAT_COPY_VALUES,&submats2[i]);CHKERRQ(ierr); 59544b287aSSatish Balay } 60544b287aSSatish Balay 61544b287aSSatish Balay /* Create submats again */ 62544b287aSSatish Balay ierr = MatCreateSubMatrices(A,n,irow,icol,MAT_REUSE_MATRIX,&submats);CHKERRQ(ierr); 63544b287aSSatish Balay 64544b287aSSatish Balay /* Compare submats and submats2 */ 65544b287aSSatish Balay for (i=0; i<n; i++) { 66544b287aSSatish Balay ierr = MatEqual(submats[i],submats2[i],&same);CHKERRQ(ierr); 67544b287aSSatish Balay if (!same) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"submatrix %d is not same\n",i); 68544b287aSSatish Balay } 69544b287aSSatish Balay 70544b287aSSatish Balay ierr = MatDestroy(&A);CHKERRQ(ierr); 71544b287aSSatish Balay for (i=0; i<n; i++) { 72544b287aSSatish Balay ierr = ISDestroy(&irow[i]);CHKERRQ(ierr); 73544b287aSSatish Balay ierr = ISDestroy(&icol[i]);CHKERRQ(ierr); 74544b287aSSatish Balay } 75544b287aSSatish Balay ierr = MatDestroySubMatrices(n,&submats);CHKERRQ(ierr); 76544b287aSSatish Balay ierr = MatDestroyMatrices(n,&submats2);CHKERRQ(ierr); 77544b287aSSatish Balay ierr = PetscFree2(irow,icol);CHKERRQ(ierr); 78544b287aSSatish Balay ierr = PetscFinalize(); 79544b287aSSatish Balay return ierr; 80544b287aSSatish Balay } 81544b287aSSatish Balay 82544b287aSSatish Balay /*TEST 83544b287aSSatish Balay 84544b287aSSatish Balay test: 85544b287aSSatish Balay suffix: 1 86544b287aSSatish Balay nsize: 2 87*dfd57a17SPierre Jolivet requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) 88544b287aSSatish Balay args: -mat_type {{aij baij}} -A ${DATAFILESPATH}/matrices/CreateSubMatrices/A -row ${DATAFILESPATH}/matrices/CreateSubMatrices/row -col ${DATAFILESPATH}/matrices/CreateSubMatrices/col 89544b287aSSatish Balay 90544b287aSSatish Balay TEST*/ 91