xref: /petsc/src/mat/tests/ex249.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
185f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
19544b287aSSatish Balay 
205f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetString(NULL,NULL,"-A",matfile,sizeof(matfile),NULL));
215f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetString(NULL,NULL,"-row",rowfile,sizeof(rowfile),NULL));
225f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetString(NULL,NULL,"-col",colfile,sizeof(colfile),NULL));
23544b287aSSatish Balay 
24544b287aSSatish Balay   /* Each rank has its own files for row/col ISes */
255f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSNPrintf(rankstr,16,"-%d",rank));
265f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrlcat(rowfile,rankstr,PETSC_MAX_PATH_LEN));
275f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrlcat(colfile,rankstr,PETSC_MAX_PATH_LEN));
28544b287aSSatish Balay 
295f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,matfile,FILE_MODE_READ,&matfd));
305f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_SELF,rowfile,FILE_MODE_READ,&rowfd));
315f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_SELF,colfile,FILE_MODE_READ,&colfd));
32544b287aSSatish Balay 
335f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
345f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(A));
355f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLoad(A,matfd));
36544b287aSSatish Balay 
37544b287aSSatish Balay   /* We stored the number of ISes at the beginning of rowfd */
385f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryRead(rowfd,&n,1,NULL,PETSC_INT));
395f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc2(n,&irow,n,&icol));
40544b287aSSatish Balay   for (i=0; i<n; i++) {
415f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreate(PETSC_COMM_SELF,&irow[i]));
425f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreate(PETSC_COMM_SELF,&icol[i]));
435f80ce2aSJacob Faibussowitsch     CHKERRQ(ISLoad(irow[i],rowfd));
445f80ce2aSJacob Faibussowitsch     CHKERRQ(ISLoad(icol[i],colfd));
45544b287aSSatish Balay   }
46544b287aSSatish Balay 
475f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&matfd));
485f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&rowfd));
495f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&colfd));
50544b287aSSatish Balay 
51544b287aSSatish Balay   /* Create submats for the first time */
525f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSubMatrices(A,n,irow,icol,MAT_INITIAL_MATRIX,&submats));
53544b287aSSatish Balay 
54544b287aSSatish Balay   /* Dup submats to submats2 for later comparison */
555f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(n,&submats2));
56544b287aSSatish Balay   for (i=0; i<n; i++) {
575f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDuplicate(submats[i],MAT_COPY_VALUES,&submats2[i]));
58544b287aSSatish Balay   }
59544b287aSSatish Balay 
60544b287aSSatish Balay   /* Create submats again */
615f80ce2aSJacob Faibussowitsch   CHKERRQ(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++) {
655f80ce2aSJacob Faibussowitsch     CHKERRQ(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 
695f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
70544b287aSSatish Balay   for (i=0; i<n; i++) {
715f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&irow[i]));
725f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&icol[i]));
73544b287aSSatish Balay   }
745f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroySubMatrices(n,&submats));
755f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroyMatrices(n,&submats2));
765f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree2(irow,icol));
77*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
78*b122ec5aSJacob 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