1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests MatMPIBAIJSetPreallocationCSR()\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown /*T 5c4762a1bSJed Brown Concepts: partitioning 6c4762a1bSJed Brown Processors: 4 7c4762a1bSJed Brown T*/ 8c4762a1bSJed Brown 9c4762a1bSJed Brown /* 10c4762a1bSJed Brown Include "petscmat.h" so that we can use matrices. Note that this file 11c4762a1bSJed Brown automatically includes: 12c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 13c4762a1bSJed Brown petscmat.h - matrices 14c4762a1bSJed Brown petscis.h - index sets 15c4762a1bSJed Brown petscviewer.h - viewers 16c4762a1bSJed Brown */ 17c4762a1bSJed Brown #include <petscmat.h> 18c4762a1bSJed Brown 19c4762a1bSJed Brown int main(int argc,char **args) 20c4762a1bSJed Brown { 21c4762a1bSJed Brown Mat A; 22c4762a1bSJed Brown PetscInt *ia,*ja, bs = 2; 23c4762a1bSJed Brown PetscInt N = 9, n; 24c4762a1bSJed Brown PetscInt rstart, rend, row, col; 25c4762a1bSJed Brown PetscInt i; 26c4762a1bSJed Brown PetscMPIInt rank,size; 27c4762a1bSJed Brown Vec v; 28c4762a1bSJed Brown 299566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 309566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 31*be096a46SBarry Smith PetscCheck(size < 5,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Can only use at most 4 processors."); 329566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 33c4762a1bSJed Brown 34c4762a1bSJed Brown /* Get a partition range based on the vector size */ 359566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(PETSC_COMM_WORLD, PETSC_DECIDE, N, &v)); 369566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(v, &n)); 379566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(v, &rstart, &rend)); 389566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v)); 39c4762a1bSJed Brown 409566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n+1,&ia)); 419566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(3*n,&ja)); 42c4762a1bSJed Brown 43c4762a1bSJed Brown /* Construct a tri-diagonal CSR indexing */ 44c4762a1bSJed Brown i = 1; 45c4762a1bSJed Brown ia[0] = 0; 46c4762a1bSJed Brown for (row = rstart; row < rend; row++) 47c4762a1bSJed Brown { 48c4762a1bSJed Brown ia[i] = ia[i-1]; 49c4762a1bSJed Brown 50c4762a1bSJed Brown /* diagonal */ 51c4762a1bSJed Brown col = row; 52c4762a1bSJed Brown { 53c4762a1bSJed Brown ja[ia[i]] = col; 54c4762a1bSJed Brown ia[i]++; 55c4762a1bSJed Brown } 56c4762a1bSJed Brown 57c4762a1bSJed Brown /* lower diagonal */ 58c4762a1bSJed Brown col = row-1; 59c4762a1bSJed Brown if (col >= 0) 60c4762a1bSJed Brown { 61c4762a1bSJed Brown ja[ia[i]] = col; 62c4762a1bSJed Brown ia[i]++; 63c4762a1bSJed Brown } 64c4762a1bSJed Brown 65c4762a1bSJed Brown /* upper diagonal */ 66c4762a1bSJed Brown col = row+1; 67c4762a1bSJed Brown if (col < N) 68c4762a1bSJed Brown { 69c4762a1bSJed Brown ja[ia[i]] = col; 70c4762a1bSJed Brown ia[i]++; 71c4762a1bSJed Brown } 72c4762a1bSJed Brown i++; 73c4762a1bSJed Brown } 74c4762a1bSJed Brown 759566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 769566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE)); 779566063dSJacob Faibussowitsch PetscCall(MatSetType(A,MATMPIAIJ)); 789566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocationCSR(A, ia, ja, NULL)); 799566063dSJacob Faibussowitsch PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 809566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 81c4762a1bSJed Brown 829566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 839566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, bs*n, bs*n, PETSC_DETERMINE, PETSC_DETERMINE)); 849566063dSJacob Faibussowitsch PetscCall(MatSetType(A,MATMPIBAIJ)); 859566063dSJacob Faibussowitsch PetscCall(MatMPIBAIJSetPreallocationCSR(A, bs, ia, ja, NULL)); 869566063dSJacob Faibussowitsch PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 879566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 88c4762a1bSJed Brown 899566063dSJacob Faibussowitsch PetscCall(PetscFree(ia)); 909566063dSJacob Faibussowitsch PetscCall(PetscFree(ja)); 919566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 92b122ec5aSJacob Faibussowitsch return 0; 93c4762a1bSJed Brown } 94c4762a1bSJed Brown 95c4762a1bSJed Brown /*TEST 96c4762a1bSJed Brown 97c4762a1bSJed Brown test: 98c4762a1bSJed Brown nsize: 4 99c4762a1bSJed Brown 100c4762a1bSJed Brown TEST*/ 101