1 2 static char help[] = "Tests MatMPIBAIJSetPreallocationCSR()\n\n"; 3 4 /*T 5 Concepts: partitioning 6 Processors: 4 7 T*/ 8 9 /* 10 Include "petscmat.h" so that we can use matrices. Note that this file 11 automatically includes: 12 petscsys.h - base PETSc routines petscvec.h - vectors 13 petscmat.h - matrices 14 petscis.h - index sets 15 petscviewer.h - viewers 16 */ 17 #include <petscmat.h> 18 19 int main(int argc,char **args) 20 { 21 Mat A; 22 PetscInt *ia,*ja, bs = 2; 23 PetscInt N = 9, n; 24 PetscInt rstart, rend, row, col; 25 PetscInt i; 26 PetscErrorCode ierr; 27 PetscMPIInt rank,size; 28 Vec v; 29 30 ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 31 CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 32 PetscCheckFalse(size > 4,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Can only use at most 4 processors."); 33 CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 34 35 /* Get a partition range based on the vector size */ 36 CHKERRQ(VecCreateMPI(PETSC_COMM_WORLD, PETSC_DECIDE, N, &v)); 37 CHKERRQ(VecGetLocalSize(v, &n)); 38 CHKERRQ(VecGetOwnershipRange(v, &rstart, &rend)); 39 CHKERRQ(VecDestroy(&v)); 40 41 CHKERRQ(PetscMalloc1(n+1,&ia)); 42 CHKERRQ(PetscMalloc1(3*n,&ja)); 43 44 /* Construct a tri-diagonal CSR indexing */ 45 i = 1; 46 ia[0] = 0; 47 for (row = rstart; row < rend; row++) 48 { 49 ia[i] = ia[i-1]; 50 51 /* diagonal */ 52 col = row; 53 { 54 ja[ia[i]] = col; 55 ia[i]++; 56 } 57 58 /* lower diagonal */ 59 col = row-1; 60 if (col >= 0) 61 { 62 ja[ia[i]] = col; 63 ia[i]++; 64 } 65 66 /* upper diagonal */ 67 col = row+1; 68 if (col < N) 69 { 70 ja[ia[i]] = col; 71 ia[i]++; 72 } 73 i++; 74 } 75 76 CHKERRQ(MatCreate(PETSC_COMM_WORLD, &A)); 77 CHKERRQ(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE)); 78 CHKERRQ(MatSetType(A,MATMPIAIJ)); 79 CHKERRQ(MatMPIAIJSetPreallocationCSR(A, ia, ja, NULL)); 80 CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 81 CHKERRQ(MatDestroy(&A)); 82 83 CHKERRQ(MatCreate(PETSC_COMM_WORLD, &A)); 84 CHKERRQ(MatSetSizes(A, bs*n, bs*n, PETSC_DETERMINE, PETSC_DETERMINE)); 85 CHKERRQ(MatSetType(A,MATMPIBAIJ)); 86 CHKERRQ(MatMPIBAIJSetPreallocationCSR(A, bs, ia, ja, NULL)); 87 CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 88 CHKERRQ(MatDestroy(&A)); 89 90 CHKERRQ(PetscFree(ia)); 91 CHKERRQ(PetscFree(ja)); 92 ierr = PetscFinalize(); 93 return ierr; 94 } 95 96 /*TEST 97 98 test: 99 nsize: 4 100 101 TEST*/ 102