1c4762a1bSJed Brown static char help[] = "Tests MatMPIAIJSetPreallocationCSR()\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown /* 4c4762a1bSJed Brown Include "petscmat.h" so that we can use matrices. Note that this file 5c4762a1bSJed Brown automatically includes: 6c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 7c4762a1bSJed Brown petscmat.h - matrices 8c4762a1bSJed Brown petscis.h - index sets 9c4762a1bSJed Brown petscviewer.h - viewers 10c4762a1bSJed Brown */ 11c4762a1bSJed Brown #include <petscmat.h> 12c4762a1bSJed Brown 13d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 14d71ae5a4SJacob Faibussowitsch { 15c4762a1bSJed Brown Mat A; 16c4762a1bSJed Brown PetscInt *ia, *ja; 17c4762a1bSJed Brown PetscMPIInt rank, size; 18c4762a1bSJed Brown 19327415f7SBarry Smith PetscFunctionBeginUser; 20*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help)); 219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 22be096a46SBarry Smith PetscCheck(size == 4, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Must run with 4 processors"); 239566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 24c4762a1bSJed Brown 259566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(5, &ia)); 269566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(16, &ja)); 27dd400576SPatrick Sanan if (rank == 0) { 289371c9d4SSatish Balay ja[0] = 1; 299371c9d4SSatish Balay ja[1] = 4; 309371c9d4SSatish Balay ja[2] = 0; 319371c9d4SSatish Balay ja[3] = 2; 329371c9d4SSatish Balay ja[4] = 5; 339371c9d4SSatish Balay ja[5] = 1; 349371c9d4SSatish Balay ja[6] = 3; 359371c9d4SSatish Balay ja[7] = 6; 369371c9d4SSatish Balay ja[8] = 2; 379371c9d4SSatish Balay ja[9] = 7; 389371c9d4SSatish Balay ia[0] = 0; 399371c9d4SSatish Balay ia[1] = 2; 409371c9d4SSatish Balay ia[2] = 5; 419371c9d4SSatish Balay ia[3] = 8; 429371c9d4SSatish Balay ia[4] = 10; 43c4762a1bSJed Brown } else if (rank == 1) { 449371c9d4SSatish Balay ja[0] = 0; 459371c9d4SSatish Balay ja[1] = 5; 469371c9d4SSatish Balay ja[2] = 8; 479371c9d4SSatish Balay ja[3] = 1; 489371c9d4SSatish Balay ja[4] = 4; 499371c9d4SSatish Balay ja[5] = 6; 509371c9d4SSatish Balay ja[6] = 9; 519371c9d4SSatish Balay ja[7] = 2; 529371c9d4SSatish Balay ja[8] = 5; 539371c9d4SSatish Balay ja[9] = 7; 549371c9d4SSatish Balay ja[10] = 10; 559371c9d4SSatish Balay ja[11] = 3; 569371c9d4SSatish Balay ja[12] = 6; 579371c9d4SSatish Balay ja[13] = 11; 589371c9d4SSatish Balay ia[0] = 0; 599371c9d4SSatish Balay ia[1] = 3; 609371c9d4SSatish Balay ia[2] = 7; 619371c9d4SSatish Balay ia[3] = 11; 629371c9d4SSatish Balay ia[4] = 14; 63c4762a1bSJed Brown } else if (rank == 2) { 649371c9d4SSatish Balay ja[0] = 4; 659371c9d4SSatish Balay ja[1] = 9; 669371c9d4SSatish Balay ja[2] = 12; 679371c9d4SSatish Balay ja[3] = 5; 689371c9d4SSatish Balay ja[4] = 8; 699371c9d4SSatish Balay ja[5] = 10; 709371c9d4SSatish Balay ja[6] = 13; 719371c9d4SSatish Balay ja[7] = 6; 729371c9d4SSatish Balay ja[8] = 9; 739371c9d4SSatish Balay ja[9] = 11; 749371c9d4SSatish Balay ja[10] = 14; 759371c9d4SSatish Balay ja[11] = 7; 769371c9d4SSatish Balay ja[12] = 10; 779371c9d4SSatish Balay ja[13] = 15; 789371c9d4SSatish Balay ia[0] = 0; 799371c9d4SSatish Balay ia[1] = 3; 809371c9d4SSatish Balay ia[2] = 7; 819371c9d4SSatish Balay ia[3] = 11; 829371c9d4SSatish Balay ia[4] = 14; 83c4762a1bSJed Brown } else { 849371c9d4SSatish Balay ja[0] = 8; 859371c9d4SSatish Balay ja[1] = 13; 869371c9d4SSatish Balay ja[2] = 9; 879371c9d4SSatish Balay ja[3] = 12; 889371c9d4SSatish Balay ja[4] = 14; 899371c9d4SSatish Balay ja[5] = 10; 909371c9d4SSatish Balay ja[6] = 13; 919371c9d4SSatish Balay ja[7] = 15; 929371c9d4SSatish Balay ja[8] = 11; 939371c9d4SSatish Balay ja[9] = 14; 949371c9d4SSatish Balay ia[0] = 0; 959371c9d4SSatish Balay ia[1] = 2; 969371c9d4SSatish Balay ia[2] = 5; 979371c9d4SSatish Balay ia[3] = 8; 989371c9d4SSatish Balay ia[4] = 10; 99c4762a1bSJed Brown } 100c4762a1bSJed Brown 1019566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 1029566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, 4, 4, 16, 16)); 1039566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATMPIAIJ)); 1049566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocationCSR(A, ia, ja, NULL)); 1059566063dSJacob Faibussowitsch PetscCall(PetscFree(ia)); 1069566063dSJacob Faibussowitsch PetscCall(PetscFree(ja)); 1079566063dSJacob Faibussowitsch PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); 1089566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1099566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 110b122ec5aSJacob Faibussowitsch return 0; 111c4762a1bSJed Brown } 112c4762a1bSJed Brown 113c4762a1bSJed Brown /*TEST 114c4762a1bSJed Brown 115c4762a1bSJed Brown test: 116c4762a1bSJed Brown nsize: 4 117c4762a1bSJed Brown 118c4762a1bSJed Brown TEST*/ 119