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 PetscErrorCode ierr; 27c4762a1bSJed Brown PetscMPIInt rank,size; 28c4762a1bSJed Brown Vec v; 29c4762a1bSJed Brown 30c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 31*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 322c71b3e2SJacob Faibussowitsch PetscCheckFalse(size > 4,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Can only use at most 4 processors."); 33*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 34c4762a1bSJed Brown 35c4762a1bSJed Brown /* Get a partition range based on the vector size */ 36*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateMPI(PETSC_COMM_WORLD, PETSC_DECIDE, N, &v)); 37*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(v, &n)); 38*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetOwnershipRange(v, &rstart, &rend)); 39*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&v)); 40c4762a1bSJed Brown 41*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(n+1,&ia)); 42*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(3*n,&ja)); 43c4762a1bSJed Brown 44c4762a1bSJed Brown /* Construct a tri-diagonal CSR indexing */ 45c4762a1bSJed Brown i = 1; 46c4762a1bSJed Brown ia[0] = 0; 47c4762a1bSJed Brown for (row = rstart; row < rend; row++) 48c4762a1bSJed Brown { 49c4762a1bSJed Brown ia[i] = ia[i-1]; 50c4762a1bSJed Brown 51c4762a1bSJed Brown /* diagonal */ 52c4762a1bSJed Brown col = row; 53c4762a1bSJed Brown { 54c4762a1bSJed Brown ja[ia[i]] = col; 55c4762a1bSJed Brown ia[i]++; 56c4762a1bSJed Brown } 57c4762a1bSJed Brown 58c4762a1bSJed Brown /* lower diagonal */ 59c4762a1bSJed Brown col = row-1; 60c4762a1bSJed Brown if (col >= 0) 61c4762a1bSJed Brown { 62c4762a1bSJed Brown ja[ia[i]] = col; 63c4762a1bSJed Brown ia[i]++; 64c4762a1bSJed Brown } 65c4762a1bSJed Brown 66c4762a1bSJed Brown /* upper diagonal */ 67c4762a1bSJed Brown col = row+1; 68c4762a1bSJed Brown if (col < N) 69c4762a1bSJed Brown { 70c4762a1bSJed Brown ja[ia[i]] = col; 71c4762a1bSJed Brown ia[i]++; 72c4762a1bSJed Brown } 73c4762a1bSJed Brown i++; 74c4762a1bSJed Brown } 75c4762a1bSJed Brown 76*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD, &A)); 77*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE)); 78*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A,MATMPIAIJ)); 79*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMPIAIJSetPreallocationCSR(A, ia, ja, NULL)); 80*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 81*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 82c4762a1bSJed Brown 83*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD, &A)); 84*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A, bs*n, bs*n, PETSC_DETERMINE, PETSC_DETERMINE)); 85*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A,MATMPIBAIJ)); 86*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMPIBAIJSetPreallocationCSR(A, bs, ia, ja, NULL)); 87*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 88*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 89c4762a1bSJed Brown 90*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(ia)); 91*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(ja)); 92c4762a1bSJed Brown ierr = PetscFinalize(); 93c4762a1bSJed Brown return ierr; 94c4762a1bSJed Brown } 95c4762a1bSJed Brown 96c4762a1bSJed Brown /*TEST 97c4762a1bSJed Brown 98c4762a1bSJed Brown test: 99c4762a1bSJed Brown nsize: 4 100c4762a1bSJed Brown 101c4762a1bSJed Brown TEST*/ 102