1c4762a1bSJed Brown static char help[] = "Tests MatCreateSubmatrix() in parallel."; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscmat.h> 4c4762a1bSJed Brown 5d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 6d71ae5a4SJacob Faibussowitsch { 7c4762a1bSJed Brown Mat C, A; 8c4762a1bSJed Brown PetscInt i, j, m = 3, n = 2, rstart, rend; 9c4762a1bSJed Brown PetscMPIInt size, rank; 10c4762a1bSJed Brown PetscScalar v; 11c4762a1bSJed Brown IS isrow, iscol; 1230203840SJunchao Zhang PetscBool test_matmatmult = PETSC_FALSE; 13c4762a1bSJed Brown 14327415f7SBarry Smith PetscFunctionBeginUser; 15*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help)); 1630203840SJunchao Zhang PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_matmatmult", &test_matmatmult, NULL)); 1730203840SJunchao Zhang 189566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 199566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 20c4762a1bSJed Brown n = 2 * size; 21c4762a1bSJed Brown 229566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &C)); 239566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, m * n, m * n)); 249566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C)); 259566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 26c4762a1bSJed Brown 27c4762a1bSJed Brown /* 28c4762a1bSJed Brown This is JUST to generate a nice test matrix, all processors fill up 29c4762a1bSJed Brown the entire matrix. This is not something one would ever do in practice. 30c4762a1bSJed Brown */ 319566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(C, &rstart, &rend)); 32c4762a1bSJed Brown for (i = rstart; i < rend; i++) { 33c4762a1bSJed Brown for (j = 0; j < m * n; j++) { 34c4762a1bSJed Brown v = i + j + 1; 359566063dSJacob Faibussowitsch PetscCall(MatSetValues(C, 1, &i, 1, &j, &v, INSERT_VALUES)); 36c4762a1bSJed Brown } 37c4762a1bSJed Brown } 38c4762a1bSJed Brown 399566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 409566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 419566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_COMMON)); 429566063dSJacob Faibussowitsch PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD)); 43c4762a1bSJed Brown 44c4762a1bSJed Brown /* 45c4762a1bSJed Brown Generate a new matrix consisting of every second row and column of 46c4762a1bSJed Brown the original matrix 47c4762a1bSJed Brown */ 489566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(C, &rstart, &rend)); 49c4762a1bSJed Brown /* Create parallel IS with the rows we want on THIS processor */ 509566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD, (rend - rstart) / 2, rstart, 2, &isrow)); 51c4762a1bSJed Brown /* Create parallel IS with the rows we want on THIS processor (same as rows for now) */ 529566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD, (rend - rstart) / 2, rstart, 2, &iscol)); 53c4762a1bSJed Brown 549566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(C, isrow, iscol, MAT_INITIAL_MATRIX, &A)); 559566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(C, isrow, iscol, MAT_REUSE_MATRIX, &A)); 569566063dSJacob Faibussowitsch PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); 57c4762a1bSJed Brown 5830203840SJunchao Zhang if (test_matmatmult) { 5930203840SJunchao Zhang PetscCall(MatDestroy(&C)); 60fb842aefSJose E. Roman PetscCall(MatMatMult(A, A, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &C)); 6130203840SJunchao Zhang PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD)); 6230203840SJunchao Zhang } 6330203840SJunchao Zhang 649566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrow)); 659566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscol)); 669566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 679566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 689566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 69b122ec5aSJacob Faibussowitsch return 0; 70c4762a1bSJed Brown } 71c4762a1bSJed Brown 72c4762a1bSJed Brown /*TEST 73c4762a1bSJed Brown 74c4762a1bSJed Brown test: 75c4762a1bSJed Brown 76c4762a1bSJed Brown test: 77c4762a1bSJed Brown suffix: 2 78c4762a1bSJed Brown nsize: 3 79c4762a1bSJed Brown 80c4762a1bSJed Brown test: 81c4762a1bSJed Brown suffix: 2_baij 82c4762a1bSJed Brown nsize: 3 83c4762a1bSJed Brown args: -mat_type baij 84c4762a1bSJed Brown 85c4762a1bSJed Brown test: 86c4762a1bSJed Brown suffix: 2_sbaij 87c4762a1bSJed Brown nsize: 3 88c4762a1bSJed Brown args: -mat_type sbaij 89c4762a1bSJed Brown 90c4762a1bSJed Brown test: 91c4762a1bSJed Brown suffix: baij 92c4762a1bSJed Brown args: -mat_type baij 93c4762a1bSJed Brown output_file: output/ex59_1_baij.out 94c4762a1bSJed Brown 95c4762a1bSJed Brown test: 96c4762a1bSJed Brown suffix: sbaij 97c4762a1bSJed Brown args: -mat_type sbaij 98c4762a1bSJed Brown output_file: output/ex59_1_sbaij.out 99c4762a1bSJed Brown 10030203840SJunchao Zhang test: 10130203840SJunchao Zhang suffix: kok 10230203840SJunchao Zhang nsize: 3 10330203840SJunchao Zhang requires: kokkos_kernels 10430203840SJunchao Zhang args: -mat_type aijkokkos -test_matmatmult 10530203840SJunchao Zhang filter: grep -v -i type 10630203840SJunchao Zhang output_file: output/ex59_kok.out 10730203840SJunchao Zhang 108c4762a1bSJed Brown TEST*/ 109