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 14*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,NULL,help)); 15*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 162c71b3e2SJacob Faibussowitsch PetscCheckFalse(size > 1,PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"Only coded for one process"); 17*9566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); 18*9566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,1)); 19*9566063dSJacob Faibussowitsch PetscCall(MatSetBlockSizes(A,2,1)); 20*9566063dSJacob Faibussowitsch PetscCall(MatSetType(A,MATAIJ)); 21*9566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,2,1,indices,PETSC_COPY_VALUES,&rmap)); 22*9566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,1,indices,PETSC_COPY_VALUES,&cmap)); 23*9566063dSJacob Faibussowitsch PetscCall(MatSetLocalToGlobalMapping(A,rmap,cmap)); 24*9566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(&rmap)); 25*9566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(&cmap)); 26*9566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,2,indices,PETSC_COPY_VALUES,&rmap)); 27*9566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,3,indices,PETSC_COPY_VALUES,&cmap)); 28*9566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&B)); 29*9566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,2,3)); 30*9566063dSJacob Faibussowitsch PetscCall(MatSetBlockSizes(B,1,1)); 31*9566063dSJacob Faibussowitsch PetscCall(MatSetType(B,MATAIJ)); 32*9566063dSJacob Faibussowitsch PetscCall(MatSetLocalToGlobalMapping(B,rmap,cmap)); 33*9566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(&rmap)); 34*9566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(&cmap)); 35*9566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 36*9566063dSJacob Faibussowitsch PetscCall(MatSetUp(B)); 3773b6653fSLawrence Mitchell mats[0] = A; 3873b6653fSLawrence Mitchell mats[1] = B; 3973b6653fSLawrence Mitchell mats[2] = A; 4073b6653fSLawrence Mitchell mats[3] = NULL; 4173b6653fSLawrence Mitchell mats[4] = B; 4273b6653fSLawrence Mitchell mats[5] = A; 43*9566063dSJacob Faibussowitsch PetscCall(MatCreateNest(PETSC_COMM_WORLD,2,NULL,3,NULL,mats,&C)); 44*9566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 45*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 46*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 47*9566063dSJacob Faibussowitsch PetscCall(MatView(C,NULL)); 48*9566063dSJacob Faibussowitsch PetscCall(MatNestGetISs(C,rows,NULL)); 4973b6653fSLawrence Mitchell for (i=0; i<2; i++) { 5073b6653fSLawrence Mitchell Mat submat; 5173b6653fSLawrence Mitchell IS cols[3]; 5273b6653fSLawrence Mitchell PetscInt j; 53*9566063dSJacob Faibussowitsch PetscCall(ISView(rows[i],NULL)); 54*9566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(C,rows[i],NULL,MAT_INITIAL_MATRIX,&submat)); 55*9566063dSJacob Faibussowitsch PetscCall(MatView(submat,NULL)); 56*9566063dSJacob Faibussowitsch PetscCall(MatNestGetISs(submat,NULL,cols)); 5773b6653fSLawrence Mitchell for (j=0; j<3; j++) { 58*9566063dSJacob Faibussowitsch PetscCall(ISView(cols[j],NULL)); 5973b6653fSLawrence Mitchell } 60*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&submat)); 6173b6653fSLawrence Mitchell } 62*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 63*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 64*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 65*9566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 66b122ec5aSJacob Faibussowitsch return 0; 6773b6653fSLawrence Mitchell } 6873b6653fSLawrence Mitchell 6973b6653fSLawrence Mitchell /*TEST 7073b6653fSLawrence Mitchell 7173b6653fSLawrence Mitchell test: 7273b6653fSLawrence Mitchell 7373b6653fSLawrence Mitchell TEST*/ 74