1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests late MatSetBlockSizes.\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown int main(int argc,char **args) 7c4762a1bSJed Brown { 8c4762a1bSJed Brown Mat A; 9c4762a1bSJed Brown Vec x[4]; 10c4762a1bSJed Brown IS is; 11c4762a1bSJed Brown ISLocalToGlobalMapping rmap,cmap; 12c4762a1bSJed Brown PetscInt bs[4],l2gbs[4],rbs,cbs,l2grbs,l2gcbs,i; 13c4762a1bSJed Brown 14*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 15*9566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); 16*9566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A,12,12,PETSC_DECIDE,PETSC_DECIDE)); 17*9566063dSJacob Faibussowitsch PetscCall(MatSetType(A,MATAIJ)); 18*9566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,12,0,1,&is)); 19*9566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreateIS(is,&rmap)); 20*9566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingSetBlockSize(rmap,2)); 21*9566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreateIS(is,&cmap)); 22*9566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingSetBlockSize(cmap,2)); 23c4762a1bSJed Brown 24*9566063dSJacob Faibussowitsch PetscCall(MatSetLocalToGlobalMapping(A,rmap,cmap)); 25*9566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(&rmap)); 26*9566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(&cmap)); 27*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 28*9566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 29c4762a1bSJed Brown 30*9566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A,&x[1],&x[0])); 31*9566063dSJacob Faibussowitsch PetscCall(MatSetBlockSizes(A,6,3)); 32*9566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A,&x[3],&x[2])); 33c4762a1bSJed Brown for (i=0;i<4;i++) { 34c4762a1bSJed Brown ISLocalToGlobalMapping l2g; 35c4762a1bSJed Brown 36*9566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(x[i],&bs[i])); 37*9566063dSJacob Faibussowitsch PetscCall(VecGetLocalToGlobalMapping(x[i],&l2g)); 38*9566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockSize(l2g,&l2gbs[i])); 39*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x[i])); 40c4762a1bSJed Brown } 41*9566063dSJacob Faibussowitsch PetscCall(MatGetBlockSizes(A,&rbs,&cbs)); 42*9566063dSJacob Faibussowitsch PetscCall(MatGetLocalToGlobalMapping(A,&rmap,&cmap)); 43*9566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockSize(rmap,&l2grbs)); 44*9566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetBlockSize(cmap,&l2gcbs)); 45*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 46*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Mat Block sizes: %" PetscInt_FMT " %" PetscInt_FMT " (l2g %" PetscInt_FMT " %" PetscInt_FMT ")\n",rbs,cbs,l2grbs,l2gcbs)); 47*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Vec Block sizes: %" PetscInt_FMT " %" PetscInt_FMT " (l2g %" PetscInt_FMT " %" PetscInt_FMT ")\n",bs[0],bs[1],l2gbs[0],l2gbs[1])); 48*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Vec Block sizes: %" PetscInt_FMT " %" PetscInt_FMT " (l2g %" PetscInt_FMT " %" PetscInt_FMT ")\n",bs[2],bs[3],l2gbs[2],l2gbs[3])); 49*9566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 50b122ec5aSJacob Faibussowitsch return 0; 51c4762a1bSJed Brown } 52c4762a1bSJed Brown 53c4762a1bSJed Brown /*TEST 54c4762a1bSJed Brown 55c4762a1bSJed Brown test: 56c4762a1bSJed Brown nsize: 2 57c4762a1bSJed Brown 58c4762a1bSJed Brown TEST*/ 59