1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests MatCreateSubmatrix() in parallel."; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 7d71ae5a4SJacob Faibussowitsch { 8c4762a1bSJed Brown Mat C, A; 9c4762a1bSJed Brown PetscInt i, j, m = 3, n = 2, rstart, rend; 10c4762a1bSJed Brown PetscMPIInt size, rank; 11c4762a1bSJed Brown PetscScalar v; 12c4762a1bSJed Brown IS isrow, iscol; 13*30203840SJunchao Zhang PetscBool test_matmatmult = PETSC_FALSE; 14c4762a1bSJed Brown 15327415f7SBarry Smith PetscFunctionBeginUser; 169566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 17*30203840SJunchao Zhang PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_matmatmult", &test_matmatmult, NULL)); 18*30203840SJunchao Zhang 199566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 209566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 21c4762a1bSJed Brown n = 2 * size; 22c4762a1bSJed Brown 239566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &C)); 249566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, m * n, m * n)); 259566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C)); 269566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 27c4762a1bSJed Brown 28c4762a1bSJed Brown /* 29c4762a1bSJed Brown This is JUST to generate a nice test matrix, all processors fill up 30c4762a1bSJed Brown the entire matrix. This is not something one would ever do in practice. 31c4762a1bSJed Brown */ 329566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(C, &rstart, &rend)); 33c4762a1bSJed Brown for (i = rstart; i < rend; i++) { 34c4762a1bSJed Brown for (j = 0; j < m * n; j++) { 35c4762a1bSJed Brown v = i + j + 1; 369566063dSJacob Faibussowitsch PetscCall(MatSetValues(C, 1, &i, 1, &j, &v, INSERT_VALUES)); 37c4762a1bSJed Brown } 38c4762a1bSJed Brown } 39c4762a1bSJed Brown 409566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 419566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 429566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_COMMON)); 439566063dSJacob Faibussowitsch PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD)); 44c4762a1bSJed Brown 45c4762a1bSJed Brown /* 46c4762a1bSJed Brown Generate a new matrix consisting of every second row and column of 47c4762a1bSJed Brown the original matrix 48c4762a1bSJed Brown */ 499566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(C, &rstart, &rend)); 50c4762a1bSJed Brown /* Create parallel IS with the rows we want on THIS processor */ 519566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD, (rend - rstart) / 2, rstart, 2, &isrow)); 52c4762a1bSJed Brown /* Create parallel IS with the rows we want on THIS processor (same as rows for now) */ 539566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD, (rend - rstart) / 2, rstart, 2, &iscol)); 54c4762a1bSJed Brown 559566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(C, isrow, iscol, MAT_INITIAL_MATRIX, &A)); 569566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(C, isrow, iscol, MAT_REUSE_MATRIX, &A)); 579566063dSJacob Faibussowitsch PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); 58c4762a1bSJed Brown 59*30203840SJunchao Zhang if (test_matmatmult) { 60*30203840SJunchao Zhang PetscCall(MatDestroy(&C)); 61*30203840SJunchao Zhang PetscCall(MatMatMult(A, A, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &C)); 62*30203840SJunchao Zhang PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD)); 63*30203840SJunchao Zhang } 64*30203840SJunchao Zhang 659566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrow)); 669566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscol)); 679566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 689566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 699566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 70b122ec5aSJacob Faibussowitsch return 0; 71c4762a1bSJed Brown } 72c4762a1bSJed Brown 73c4762a1bSJed Brown /*TEST 74c4762a1bSJed Brown 75c4762a1bSJed Brown test: 76c4762a1bSJed Brown 77c4762a1bSJed Brown test: 78c4762a1bSJed Brown suffix: 2 79c4762a1bSJed Brown nsize: 3 80c4762a1bSJed Brown 81c4762a1bSJed Brown test: 82c4762a1bSJed Brown suffix: 2_baij 83c4762a1bSJed Brown nsize: 3 84c4762a1bSJed Brown args: -mat_type baij 85c4762a1bSJed Brown 86c4762a1bSJed Brown test: 87c4762a1bSJed Brown suffix: 2_sbaij 88c4762a1bSJed Brown nsize: 3 89c4762a1bSJed Brown args: -mat_type sbaij 90c4762a1bSJed Brown 91c4762a1bSJed Brown test: 92c4762a1bSJed Brown suffix: baij 93c4762a1bSJed Brown args: -mat_type baij 94c4762a1bSJed Brown output_file: output/ex59_1_baij.out 95c4762a1bSJed Brown 96c4762a1bSJed Brown test: 97c4762a1bSJed Brown suffix: sbaij 98c4762a1bSJed Brown args: -mat_type sbaij 99c4762a1bSJed Brown output_file: output/ex59_1_sbaij.out 100c4762a1bSJed Brown 101*30203840SJunchao Zhang test: 102*30203840SJunchao Zhang suffix: kok 103*30203840SJunchao Zhang nsize: 3 104*30203840SJunchao Zhang requires: kokkos_kernels 105*30203840SJunchao Zhang args: -mat_type aijkokkos -test_matmatmult 106*30203840SJunchao Zhang filter: grep -v -i type 107*30203840SJunchao Zhang output_file: output/ex59_kok.out 108*30203840SJunchao Zhang 109c4762a1bSJed Brown TEST*/ 110