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