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*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 155f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 165f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A,12,12,PETSC_DECIDE,PETSC_DECIDE)); 175f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A,MATAIJ)); 185f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateStride(PETSC_COMM_WORLD,12,0,1,&is)); 195f80ce2aSJacob Faibussowitsch CHKERRQ(ISLocalToGlobalMappingCreateIS(is,&rmap)); 205f80ce2aSJacob Faibussowitsch CHKERRQ(ISLocalToGlobalMappingSetBlockSize(rmap,2)); 215f80ce2aSJacob Faibussowitsch CHKERRQ(ISLocalToGlobalMappingCreateIS(is,&cmap)); 225f80ce2aSJacob Faibussowitsch CHKERRQ(ISLocalToGlobalMappingSetBlockSize(cmap,2)); 23c4762a1bSJed Brown 245f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetLocalToGlobalMapping(A,rmap,cmap)); 255f80ce2aSJacob Faibussowitsch CHKERRQ(ISLocalToGlobalMappingDestroy(&rmap)); 265f80ce2aSJacob Faibussowitsch CHKERRQ(ISLocalToGlobalMappingDestroy(&cmap)); 275f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&is)); 285f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(A)); 29c4762a1bSJed Brown 305f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(A,&x[1],&x[0])); 315f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetBlockSizes(A,6,3)); 325f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(A,&x[3],&x[2])); 33c4762a1bSJed Brown for (i=0;i<4;i++) { 34c4762a1bSJed Brown ISLocalToGlobalMapping l2g; 35c4762a1bSJed Brown 365f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetBlockSize(x[i],&bs[i])); 375f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalToGlobalMapping(x[i],&l2g)); 385f80ce2aSJacob Faibussowitsch CHKERRQ(ISLocalToGlobalMappingGetBlockSize(l2g,&l2gbs[i])); 395f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x[i])); 40c4762a1bSJed Brown } 415f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetBlockSizes(A,&rbs,&cbs)); 425f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalToGlobalMapping(A,&rmap,&cmap)); 435f80ce2aSJacob Faibussowitsch CHKERRQ(ISLocalToGlobalMappingGetBlockSize(rmap,&l2grbs)); 445f80ce2aSJacob Faibussowitsch CHKERRQ(ISLocalToGlobalMappingGetBlockSize(cmap,&l2gcbs)); 455f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Mat Block sizes: %" PetscInt_FMT " %" PetscInt_FMT " (l2g %" PetscInt_FMT " %" PetscInt_FMT ")\n",rbs,cbs,l2grbs,l2gcbs)); 475f80ce2aSJacob Faibussowitsch CHKERRQ(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])); 485f80ce2aSJacob Faibussowitsch CHKERRQ(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*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 50*b122ec5aSJacob 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