xref: /petsc/src/mat/tests/ex212.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1 static char help[] = "Test MatCreateSubMatrix with -mat_type nest and block sizes.\n";
2 
3 #include <petscmat.h>
4 
5 int main(int argc, char **argv)
6 {
7   Mat                    A, B, C, mats[6];
8   IS                     rows[2];
9   ISLocalToGlobalMapping cmap, rmap;
10   const PetscInt         indices[3] = {0, 1, 2};
11   PetscInt               i;
12   PetscMPIInt            size;
13 
14   CHKERRQ(PetscInitialize(&argc,&argv,NULL,help));
15   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
16   PetscCheckFalse(size > 1,PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"Only coded for one process");
17   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
18   CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,1));
19   CHKERRQ(MatSetBlockSizes(A,2,1));
20   CHKERRQ(MatSetType(A,MATAIJ));
21   CHKERRQ(ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,2,1,indices,PETSC_COPY_VALUES,&rmap));
22   CHKERRQ(ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,1,indices,PETSC_COPY_VALUES,&cmap));
23   CHKERRQ(MatSetLocalToGlobalMapping(A,rmap,cmap));
24   CHKERRQ(ISLocalToGlobalMappingDestroy(&rmap));
25   CHKERRQ(ISLocalToGlobalMappingDestroy(&cmap));
26   CHKERRQ(ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,2,indices,PETSC_COPY_VALUES,&rmap));
27   CHKERRQ(ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,3,indices,PETSC_COPY_VALUES,&cmap));
28   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&B));
29   CHKERRQ(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,2,3));
30   CHKERRQ(MatSetBlockSizes(B,1,1));
31   CHKERRQ(MatSetType(B,MATAIJ));
32   CHKERRQ(MatSetLocalToGlobalMapping(B,rmap,cmap));
33   CHKERRQ(ISLocalToGlobalMappingDestroy(&rmap));
34   CHKERRQ(ISLocalToGlobalMappingDestroy(&cmap));
35   CHKERRQ(MatSetUp(A));
36   CHKERRQ(MatSetUp(B));
37   mats[0] = A;
38   mats[1] = B;
39   mats[2] = A;
40   mats[3] = NULL;
41   mats[4] = B;
42   mats[5] = A;
43   CHKERRQ(MatCreateNest(PETSC_COMM_WORLD,2,NULL,3,NULL,mats,&C));
44   CHKERRQ(MatSetUp(C));
45   CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
46   CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
47   CHKERRQ(MatView(C,NULL));
48   CHKERRQ(MatNestGetISs(C,rows,NULL));
49   for (i=0; i<2; i++) {
50     Mat      submat;
51     IS       cols[3];
52     PetscInt j;
53     CHKERRQ(ISView(rows[i],NULL));
54     CHKERRQ(MatCreateSubMatrix(C,rows[i],NULL,MAT_INITIAL_MATRIX,&submat));
55     CHKERRQ(MatView(submat,NULL));
56     CHKERRQ(MatNestGetISs(submat,NULL,cols));
57     for (j=0; j<3; j++) {
58       CHKERRQ(ISView(cols[j],NULL));
59     }
60     CHKERRQ(MatDestroy(&submat));
61   }
62   CHKERRQ(MatDestroy(&C));
63   CHKERRQ(MatDestroy(&B));
64   CHKERRQ(MatDestroy(&A));
65   CHKERRQ(PetscFinalize());
66   return 0;
67 }
68 
69 /*TEST
70 
71     test:
72 
73 TEST*/
74