1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests MatMPIAIJSetPreallocationCSR()\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; 23c4762a1bSJed Brown PetscMPIInt rank,size; 24c4762a1bSJed Brown 25*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 265f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 272c71b3e2SJacob Faibussowitsch PetscCheckFalse(size != 4,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Must run with 4 processors"); 285f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 29c4762a1bSJed Brown 305f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(5,&ia)); 315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(16,&ja)); 32dd400576SPatrick Sanan if (rank == 0) { 33c4762a1bSJed Brown ja[0] = 1; ja[1] = 4; ja[2] = 0; ja[3] = 2; ja[4] = 5; ja[5] = 1; ja[6] = 3; ja[7] = 6; 34c4762a1bSJed Brown ja[8] = 2; ja[9] = 7; 35c4762a1bSJed Brown ia[0] = 0; ia[1] = 2; ia[2] = 5; ia[3] = 8; ia[4] = 10; 36c4762a1bSJed Brown } else if (rank == 1) { 37c4762a1bSJed Brown ja[0] = 0; ja[1] = 5; ja[2] = 8; ja[3] = 1; ja[4] = 4; ja[5] = 6; ja[6] = 9; ja[7] = 2; 38c4762a1bSJed Brown ja[8] = 5; ja[9] = 7; ja[10] = 10; ja[11] = 3; ja[12] = 6; ja[13] = 11; 39c4762a1bSJed Brown ia[0] = 0; ia[1] = 3; ia[2] = 7; ia[3] = 11; ia[4] = 14; 40c4762a1bSJed Brown } else if (rank == 2) { 41c4762a1bSJed Brown ja[0] = 4; ja[1] = 9; ja[2] = 12; ja[3] = 5; ja[4] = 8; ja[5] = 10; ja[6] = 13; ja[7] = 6; 42c4762a1bSJed Brown ja[8] = 9; ja[9] = 11; ja[10] = 14; ja[11] = 7; ja[12] = 10; ja[13] = 15; 43c4762a1bSJed Brown ia[0] = 0; ia[1] = 3; ia[2] = 7; ia[3] = 11; ia[4] = 14; 44c4762a1bSJed Brown } else { 45c4762a1bSJed Brown ja[0] = 8; ja[1] = 13; ja[2] = 9; ja[3] = 12; ja[4] = 14; ja[5] = 10; ja[6] = 13; ja[7] = 15; 46c4762a1bSJed Brown ja[8] = 11; ja[9] = 14; 47c4762a1bSJed Brown ia[0] = 0; ia[1] = 2; ia[2] = 5; ia[3] = 8; ia[4] = 10; 48c4762a1bSJed Brown } 49c4762a1bSJed Brown 505f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 515f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A,4,4,16,16)); 525f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A,MATMPIAIJ)); 535f80ce2aSJacob Faibussowitsch CHKERRQ(MatMPIAIJSetPreallocationCSR(A,ia,ja,NULL)); 545f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(ia)); 555f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(ja)); 565f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 575f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 58*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 59*b122ec5aSJacob Faibussowitsch return 0; 60c4762a1bSJed Brown } 61c4762a1bSJed Brown 62c4762a1bSJed Brown /*TEST 63c4762a1bSJed Brown 64c4762a1bSJed Brown test: 65c4762a1bSJed Brown nsize: 4 66c4762a1bSJed Brown 67c4762a1bSJed Brown TEST*/ 68