1 2 static char help[] = "Tests MatMPIAIJSetPreallocationCSR()\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; 23 PetscMPIInt rank,size; 24 25 PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 26 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 27 PetscCheckFalse(size != 4,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Must run with 4 processors"); 28 PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 29 30 PetscCall(PetscMalloc1(5,&ia)); 31 PetscCall(PetscMalloc1(16,&ja)); 32 if (rank == 0) { 33 ja[0] = 1; ja[1] = 4; ja[2] = 0; ja[3] = 2; ja[4] = 5; ja[5] = 1; ja[6] = 3; ja[7] = 6; 34 ja[8] = 2; ja[9] = 7; 35 ia[0] = 0; ia[1] = 2; ia[2] = 5; ia[3] = 8; ia[4] = 10; 36 } else if (rank == 1) { 37 ja[0] = 0; ja[1] = 5; ja[2] = 8; ja[3] = 1; ja[4] = 4; ja[5] = 6; ja[6] = 9; ja[7] = 2; 38 ja[8] = 5; ja[9] = 7; ja[10] = 10; ja[11] = 3; ja[12] = 6; ja[13] = 11; 39 ia[0] = 0; ia[1] = 3; ia[2] = 7; ia[3] = 11; ia[4] = 14; 40 } else if (rank == 2) { 41 ja[0] = 4; ja[1] = 9; ja[2] = 12; ja[3] = 5; ja[4] = 8; ja[5] = 10; ja[6] = 13; ja[7] = 6; 42 ja[8] = 9; ja[9] = 11; ja[10] = 14; ja[11] = 7; ja[12] = 10; ja[13] = 15; 43 ia[0] = 0; ia[1] = 3; ia[2] = 7; ia[3] = 11; ia[4] = 14; 44 } else { 45 ja[0] = 8; ja[1] = 13; ja[2] = 9; ja[3] = 12; ja[4] = 14; ja[5] = 10; ja[6] = 13; ja[7] = 15; 46 ja[8] = 11; ja[9] = 14; 47 ia[0] = 0; ia[1] = 2; ia[2] = 5; ia[3] = 8; ia[4] = 10; 48 } 49 50 PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); 51 PetscCall(MatSetSizes(A,4,4,16,16)); 52 PetscCall(MatSetType(A,MATMPIAIJ)); 53 PetscCall(MatMPIAIJSetPreallocationCSR(A,ia,ja,NULL)); 54 PetscCall(PetscFree(ia)); 55 PetscCall(PetscFree(ja)); 56 PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 57 PetscCall(MatDestroy(&A)); 58 PetscCall(PetscFinalize()); 59 return 0; 60 } 61 62 /*TEST 63 64 test: 65 nsize: 4 66 67 TEST*/ 68