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