xref: /petsc/src/mat/tests/ex212.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
173b6653fSLawrence Mitchell static char help[] = "Test MatCreateSubMatrix with -mat_type nest and block sizes.\n";
273b6653fSLawrence Mitchell 
373b6653fSLawrence Mitchell #include <petscmat.h>
473b6653fSLawrence Mitchell 
573b6653fSLawrence Mitchell int main(int argc, char **argv)
673b6653fSLawrence Mitchell {
773b6653fSLawrence Mitchell   Mat                    A, B, C, mats[6];
873b6653fSLawrence Mitchell   IS                     rows[2];
973b6653fSLawrence Mitchell   ISLocalToGlobalMapping cmap, rmap;
1073b6653fSLawrence Mitchell   const PetscInt         indices[3] = {0, 1, 2};
1173b6653fSLawrence Mitchell   PetscInt               i;
1273b6653fSLawrence Mitchell   PetscMPIInt            size;
1373b6653fSLawrence Mitchell   PetscErrorCode         ierr;
1473b6653fSLawrence Mitchell 
1573b6653fSLawrence Mitchell   ierr = PetscInitialize(&argc,&argv,NULL,help); if (ierr) return ierr;
16*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
172c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size > 1,PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"Only coded for one process");
18*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
19*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,1));
20*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetBlockSizes(A,2,1));
21*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(A,MATAIJ));
22*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,2,1,indices,PETSC_COPY_VALUES,&rmap));
23*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,1,indices,PETSC_COPY_VALUES,&cmap));
24*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetLocalToGlobalMapping(A,rmap,cmap));
25*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISLocalToGlobalMappingDestroy(&rmap));
26*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISLocalToGlobalMappingDestroy(&cmap));
27*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,2,indices,PETSC_COPY_VALUES,&rmap));
28*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,3,indices,PETSC_COPY_VALUES,&cmap));
29*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&B));
30*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,2,3));
31*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetBlockSizes(B,1,1));
32*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(B,MATAIJ));
33*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetLocalToGlobalMapping(B,rmap,cmap));
34*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISLocalToGlobalMappingDestroy(&rmap));
35*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISLocalToGlobalMappingDestroy(&cmap));
36*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(A));
37*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(B));
3873b6653fSLawrence Mitchell   mats[0] = A;
3973b6653fSLawrence Mitchell   mats[1] = B;
4073b6653fSLawrence Mitchell   mats[2] = A;
4173b6653fSLawrence Mitchell   mats[3] = NULL;
4273b6653fSLawrence Mitchell   mats[4] = B;
4373b6653fSLawrence Mitchell   mats[5] = A;
44*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateNest(PETSC_COMM_WORLD,2,NULL,3,NULL,mats,&C));
45*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(C));
46*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
47*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
48*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(C,NULL));
49*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatNestGetISs(C,rows,NULL));
5073b6653fSLawrence Mitchell   for (i=0; i<2; i++) {
5173b6653fSLawrence Mitchell     Mat      submat;
5273b6653fSLawrence Mitchell     IS       cols[3];
5373b6653fSLawrence Mitchell     PetscInt j;
54*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISView(rows[i],NULL));
55*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateSubMatrix(C,rows[i],NULL,MAT_INITIAL_MATRIX,&submat));
56*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatView(submat,NULL));
57*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatNestGetISs(submat,NULL,cols));
5873b6653fSLawrence Mitchell     for (j=0; j<3; j++) {
59*5f80ce2aSJacob Faibussowitsch       CHKERRQ(ISView(cols[j],NULL));
6073b6653fSLawrence Mitchell     }
61*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&submat));
6273b6653fSLawrence Mitchell   }
63*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
64*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
65*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
6673b6653fSLawrence Mitchell   ierr = PetscFinalize();
6773b6653fSLawrence Mitchell   return ierr;
6873b6653fSLawrence Mitchell }
6973b6653fSLawrence Mitchell 
7073b6653fSLawrence Mitchell /*TEST
7173b6653fSLawrence Mitchell 
7273b6653fSLawrence Mitchell     test:
7373b6653fSLawrence Mitchell 
7473b6653fSLawrence Mitchell TEST*/
75