xref: /petsc/src/mat/tests/ex212.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
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;
16ffc4695bSBarry Smith   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
17*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size > 1,PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"Only coded for one process");
1873b6653fSLawrence Mitchell   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
1973b6653fSLawrence Mitchell   ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,1);CHKERRQ(ierr);
2073b6653fSLawrence Mitchell   ierr = MatSetBlockSizes(A,2,1);CHKERRQ(ierr);
2173b6653fSLawrence Mitchell   ierr = MatSetType(A,MATAIJ);CHKERRQ(ierr);
2273b6653fSLawrence Mitchell   ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,2,1,indices,PETSC_COPY_VALUES,&rmap);CHKERRQ(ierr);
2373b6653fSLawrence Mitchell   ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,1,indices,PETSC_COPY_VALUES,&cmap);CHKERRQ(ierr);
2473b6653fSLawrence Mitchell   ierr = MatSetLocalToGlobalMapping(A,rmap,cmap);CHKERRQ(ierr);
2573b6653fSLawrence Mitchell   ierr = ISLocalToGlobalMappingDestroy(&rmap);CHKERRQ(ierr);
2673b6653fSLawrence Mitchell   ierr = ISLocalToGlobalMappingDestroy(&cmap);CHKERRQ(ierr);
2773b6653fSLawrence Mitchell   ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,2,indices,PETSC_COPY_VALUES,&rmap);CHKERRQ(ierr);
2873b6653fSLawrence Mitchell   ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,3,indices,PETSC_COPY_VALUES,&cmap);CHKERRQ(ierr);
2973b6653fSLawrence Mitchell   ierr = MatCreate(PETSC_COMM_WORLD,&B);CHKERRQ(ierr);
3073b6653fSLawrence Mitchell   ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,2,3);CHKERRQ(ierr);
3173b6653fSLawrence Mitchell   ierr = MatSetBlockSizes(B,1,1);CHKERRQ(ierr);
3273b6653fSLawrence Mitchell   ierr = MatSetType(B,MATAIJ);CHKERRQ(ierr);
3373b6653fSLawrence Mitchell   ierr = MatSetLocalToGlobalMapping(B,rmap,cmap);CHKERRQ(ierr);
3473b6653fSLawrence Mitchell   ierr = ISLocalToGlobalMappingDestroy(&rmap);CHKERRQ(ierr);
3573b6653fSLawrence Mitchell   ierr = ISLocalToGlobalMappingDestroy(&cmap);CHKERRQ(ierr);
3673b6653fSLawrence Mitchell   ierr = MatSetUp(A);CHKERRQ(ierr);
3773b6653fSLawrence Mitchell   ierr = MatSetUp(B);CHKERRQ(ierr);
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;
4473b6653fSLawrence Mitchell   ierr = MatCreateNest(PETSC_COMM_WORLD,2,NULL,3,NULL,mats,&C);CHKERRQ(ierr);
4573b6653fSLawrence Mitchell   ierr = MatSetUp(C);CHKERRQ(ierr);
4673b6653fSLawrence Mitchell   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4773b6653fSLawrence Mitchell   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
4873b6653fSLawrence Mitchell   ierr = MatView(C,NULL);CHKERRQ(ierr);
4973b6653fSLawrence Mitchell   ierr = MatNestGetISs(C,rows,NULL);CHKERRQ(ierr);
5073b6653fSLawrence Mitchell   for (i=0; i<2; i++) {
5173b6653fSLawrence Mitchell     Mat      submat;
5273b6653fSLawrence Mitchell     IS       cols[3];
5373b6653fSLawrence Mitchell     PetscInt j;
5473b6653fSLawrence Mitchell     ierr = ISView(rows[i],NULL);CHKERRQ(ierr);
5573b6653fSLawrence Mitchell     ierr = MatCreateSubMatrix(C,rows[i],NULL,MAT_INITIAL_MATRIX,&submat);CHKERRQ(ierr);
5673b6653fSLawrence Mitchell     ierr = MatView(submat,NULL);CHKERRQ(ierr);
5773b6653fSLawrence Mitchell     ierr = MatNestGetISs(submat,NULL,cols);CHKERRQ(ierr);
5873b6653fSLawrence Mitchell     for (j=0; j<3; j++) {
5973b6653fSLawrence Mitchell       ierr = ISView(cols[j],NULL);CHKERRQ(ierr);
6073b6653fSLawrence Mitchell     }
6173b6653fSLawrence Mitchell     ierr = MatDestroy(&submat);CHKERRQ(ierr);
6273b6653fSLawrence Mitchell   }
6373b6653fSLawrence Mitchell   ierr = MatDestroy(&C);CHKERRQ(ierr);
6473b6653fSLawrence Mitchell   ierr = MatDestroy(&B);CHKERRQ(ierr);
6573b6653fSLawrence Mitchell   ierr = MatDestroy(&A);CHKERRQ(ierr);
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