xref: /petsc/src/mat/tests/ex249.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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;
19*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
20544b287aSSatish Balay 
21*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetString(NULL,NULL,"-A",matfile,sizeof(matfile),NULL));
22*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetString(NULL,NULL,"-row",rowfile,sizeof(rowfile),NULL));
23*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetString(NULL,NULL,"-col",colfile,sizeof(colfile),NULL));
24544b287aSSatish Balay 
25544b287aSSatish Balay   /* Each rank has its own files for row/col ISes */
26*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSNPrintf(rankstr,16,"-%d",rank));
27*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrlcat(rowfile,rankstr,PETSC_MAX_PATH_LEN));
28*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrlcat(colfile,rankstr,PETSC_MAX_PATH_LEN));
29544b287aSSatish Balay 
30*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,matfile,FILE_MODE_READ,&matfd));
31*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_SELF,rowfile,FILE_MODE_READ,&rowfd));
32*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_SELF,colfile,FILE_MODE_READ,&colfd));
33544b287aSSatish Balay 
34*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
35*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(A));
36*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLoad(A,matfd));
37544b287aSSatish Balay 
38544b287aSSatish Balay   /* We stored the number of ISes at the beginning of rowfd */
39*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryRead(rowfd,&n,1,NULL,PETSC_INT));
40*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc2(n,&irow,n,&icol));
41544b287aSSatish Balay   for (i=0; i<n; i++) {
42*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreate(PETSC_COMM_SELF,&irow[i]));
43*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreate(PETSC_COMM_SELF,&icol[i]));
44*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISLoad(irow[i],rowfd));
45*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISLoad(icol[i],colfd));
46544b287aSSatish Balay   }
47544b287aSSatish Balay 
48*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&matfd));
49*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&rowfd));
50*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&colfd));
51544b287aSSatish Balay 
52544b287aSSatish Balay   /* Create submats for the first time */
53*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSubMatrices(A,n,irow,icol,MAT_INITIAL_MATRIX,&submats));
54544b287aSSatish Balay 
55544b287aSSatish Balay   /* Dup submats to submats2 for later comparison */
56*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(n,&submats2));
57544b287aSSatish Balay   for (i=0; i<n; i++) {
58*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDuplicate(submats[i],MAT_COPY_VALUES,&submats2[i]));
59544b287aSSatish Balay   }
60544b287aSSatish Balay 
61544b287aSSatish Balay   /* Create submats again */
62*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSubMatrices(A,n,irow,icol,MAT_REUSE_MATRIX,&submats));
63544b287aSSatish Balay 
64544b287aSSatish Balay   /* Compare submats and submats2 */
65544b287aSSatish Balay   for (i=0; i<n; i++) {
66*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatEqual(submats[i],submats2[i],&same));
672c71b3e2SJacob Faibussowitsch     PetscCheckFalse(!same,PETSC_COMM_SELF,PETSC_ERR_PLIB,"submatrix %" PetscInt_FMT " is not same",i);
68544b287aSSatish Balay   }
69544b287aSSatish Balay 
70*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
71544b287aSSatish Balay   for (i=0; i<n; i++) {
72*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&irow[i]));
73*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&icol[i]));
74544b287aSSatish Balay   }
75*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroySubMatrices(n,&submats));
76*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroyMatrices(n,&submats2));
77*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree2(irow,icol));
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
87dfd57a17SPierre 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