1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests MPI parallel matrix creation. Test MatCreateRedundantMatrix() \n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6*d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 7*d71ae5a4SJacob Faibussowitsch { 8c4762a1bSJed Brown Mat C, Credundant; 9c4762a1bSJed Brown MatInfo info; 10c4762a1bSJed Brown PetscMPIInt rank, size, subsize; 11c4762a1bSJed Brown PetscInt i, j, m = 3, n = 2, low, high, iglobal; 12c4762a1bSJed Brown PetscInt Ii, J, ldim, nsubcomms; 13c4762a1bSJed Brown PetscBool flg_info, flg_mat; 14c4762a1bSJed Brown PetscScalar v, one = 1.0; 15c4762a1bSJed Brown Vec x, y; 16c4762a1bSJed Brown MPI_Comm subcomm; 17c4762a1bSJed Brown 18327415f7SBarry Smith PetscFunctionBeginUser; 199566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 209566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL)); 219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 23c4762a1bSJed Brown n = 2 * size; 24c4762a1bSJed Brown 259566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &C)); 269566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, m * n, m * n)); 279566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C)); 289566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 29c4762a1bSJed Brown 30c4762a1bSJed Brown /* Create the matrix for the five point stencil, YET AGAIN */ 31c4762a1bSJed Brown for (i = 0; i < m; i++) { 32c4762a1bSJed Brown for (j = 2 * rank; j < 2 * rank + 2; j++) { 339371c9d4SSatish Balay v = -1.0; 349371c9d4SSatish Balay Ii = j + n * i; 359371c9d4SSatish Balay if (i > 0) { 369371c9d4SSatish Balay J = Ii - n; 379371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 389371c9d4SSatish Balay } 399371c9d4SSatish Balay if (i < m - 1) { 409371c9d4SSatish Balay J = Ii + n; 419371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 429371c9d4SSatish Balay } 439371c9d4SSatish Balay if (j > 0) { 449371c9d4SSatish Balay J = Ii - 1; 459371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 469371c9d4SSatish Balay } 479371c9d4SSatish Balay if (j < n - 1) { 489371c9d4SSatish Balay J = Ii + 1; 499371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 509371c9d4SSatish Balay } 519371c9d4SSatish Balay v = 4.0; 529371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &Ii, &v, INSERT_VALUES)); 53c4762a1bSJed Brown } 54c4762a1bSJed Brown } 55c4762a1bSJed Brown 56c4762a1bSJed Brown /* Add extra elements (to illustrate variants of MatGetInfo) */ 579371c9d4SSatish Balay Ii = n; 589371c9d4SSatish Balay J = n - 2; 599371c9d4SSatish Balay v = 100.0; 609566063dSJacob Faibussowitsch PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 619371c9d4SSatish Balay Ii = n - 2; 629371c9d4SSatish Balay J = n; 639371c9d4SSatish Balay v = 100.0; 649566063dSJacob Faibussowitsch PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES)); 65c4762a1bSJed Brown 669566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 679566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 68c4762a1bSJed Brown 69c4762a1bSJed Brown /* Form vectors */ 709566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(C, &x, &y)); 719566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(x, &ldim)); 729566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(x, &low, &high)); 73c4762a1bSJed Brown for (i = 0; i < ldim; i++) { 74c4762a1bSJed Brown iglobal = i + low; 75c4762a1bSJed Brown v = one * ((PetscReal)i) + 100.0 * rank; 769566063dSJacob Faibussowitsch PetscCall(VecSetValues(x, 1, &iglobal, &v, INSERT_VALUES)); 77c4762a1bSJed Brown } 789566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(x)); 799566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(x)); 80c4762a1bSJed Brown 819566063dSJacob Faibussowitsch PetscCall(MatMult(C, x, y)); 82c4762a1bSJed Brown 839566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-view_info", &flg_info)); 84c4762a1bSJed Brown if (flg_info) { 859566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_INFO)); 869566063dSJacob Faibussowitsch PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD)); 87c4762a1bSJed Brown 889566063dSJacob Faibussowitsch PetscCall(MatGetInfo(C, MAT_GLOBAL_SUM, &info)); 899566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "matrix information (global sums):\nnonzeros = %" PetscInt_FMT ", allocated nonzeros = %" PetscInt_FMT "\n", (PetscInt)info.nz_used, (PetscInt)info.nz_allocated)); 909566063dSJacob Faibussowitsch PetscCall(MatGetInfo(C, MAT_GLOBAL_MAX, &info)); 919566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "matrix information (global max):\nnonzeros = %" PetscInt_FMT ", allocated nonzeros = %" PetscInt_FMT "\n", (PetscInt)info.nz_used, (PetscInt)info.nz_allocated)); 92c4762a1bSJed Brown } 93c4762a1bSJed Brown 949566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-view_mat", &flg_mat)); 951baa6e33SBarry Smith if (flg_mat) PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD)); 96c4762a1bSJed Brown 97c4762a1bSJed Brown /* Test MatCreateRedundantMatrix() */ 98c4762a1bSJed Brown nsubcomms = size; 999566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-nsubcomms", &nsubcomms, NULL)); 1009566063dSJacob Faibussowitsch PetscCall(MatCreateRedundantMatrix(C, nsubcomms, MPI_COMM_NULL, MAT_INITIAL_MATRIX, &Credundant)); 1019566063dSJacob Faibussowitsch PetscCall(MatCreateRedundantMatrix(C, nsubcomms, MPI_COMM_NULL, MAT_REUSE_MATRIX, &Credundant)); 102c4762a1bSJed Brown 1039566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)Credundant, &subcomm)); 1049566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(subcomm, &subsize)); 105c4762a1bSJed Brown 106c4762a1bSJed Brown if (subsize == 2 && flg_mat) { 1079566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_(subcomm), "\n[%d] Credundant:\n", rank)); 1089566063dSJacob Faibussowitsch PetscCall(MatView(Credundant, PETSC_VIEWER_STDOUT_(subcomm))); 109c4762a1bSJed Brown } 1109566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Credundant)); 111c4762a1bSJed Brown 112c4762a1bSJed Brown /* Test MatCreateRedundantMatrix() with user-provided subcomm */ 113c4762a1bSJed Brown { 114c4762a1bSJed Brown PetscSubcomm psubcomm; 115c4762a1bSJed Brown 1169566063dSJacob Faibussowitsch PetscCall(PetscSubcommCreate(PETSC_COMM_WORLD, &psubcomm)); 1179566063dSJacob Faibussowitsch PetscCall(PetscSubcommSetNumber(psubcomm, nsubcomms)); 1189566063dSJacob Faibussowitsch PetscCall(PetscSubcommSetType(psubcomm, PETSC_SUBCOMM_CONTIGUOUS)); 119c4762a1bSJed Brown /* enable runtime switch of psubcomm type, e.g., '-psubcomm_type interlaced */ 1209566063dSJacob Faibussowitsch PetscCall(PetscSubcommSetFromOptions(psubcomm)); 121c4762a1bSJed Brown 1229566063dSJacob Faibussowitsch PetscCall(MatCreateRedundantMatrix(C, nsubcomms, PetscSubcommChild(psubcomm), MAT_INITIAL_MATRIX, &Credundant)); 1239566063dSJacob Faibussowitsch PetscCall(MatCreateRedundantMatrix(C, nsubcomms, PetscSubcommChild(psubcomm), MAT_REUSE_MATRIX, &Credundant)); 124c4762a1bSJed Brown 1259566063dSJacob Faibussowitsch PetscCall(PetscSubcommDestroy(&psubcomm)); 1269566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Credundant)); 127c4762a1bSJed Brown } 128c4762a1bSJed Brown 1299566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1309566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 1319566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 1329566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 133b122ec5aSJacob Faibussowitsch return 0; 134c4762a1bSJed Brown } 135c4762a1bSJed Brown 136c4762a1bSJed Brown /*TEST 137c4762a1bSJed Brown 138c4762a1bSJed Brown test: 139c4762a1bSJed Brown nsize: 3 140c4762a1bSJed Brown args: -view_info 141c4762a1bSJed Brown 142c4762a1bSJed Brown test: 143c4762a1bSJed Brown suffix: 2 144c4762a1bSJed Brown nsize: 3 145c4762a1bSJed Brown args: -nsubcomms 2 -view_mat -psubcomm_type interlaced 146c4762a1bSJed Brown 147c4762a1bSJed Brown test: 148c4762a1bSJed Brown suffix: 3 149c4762a1bSJed Brown nsize: 3 150c4762a1bSJed Brown args: -nsubcomms 2 -view_mat -psubcomm_type contiguous 151c4762a1bSJed Brown 152c4762a1bSJed Brown test: 153c4762a1bSJed Brown suffix: 3_baij 154c4762a1bSJed Brown nsize: 3 155c4762a1bSJed Brown args: -mat_type baij -nsubcomms 2 -view_mat 156c4762a1bSJed Brown 157c4762a1bSJed Brown test: 158c4762a1bSJed Brown suffix: 3_sbaij 159c4762a1bSJed Brown nsize: 3 160c4762a1bSJed Brown args: -mat_type sbaij -nsubcomms 2 -view_mat 161c4762a1bSJed Brown 162c4762a1bSJed Brown test: 163c4762a1bSJed Brown suffix: 3_dense 164c4762a1bSJed Brown nsize: 3 165c4762a1bSJed Brown args: -mat_type dense -nsubcomms 2 -view_mat 166c4762a1bSJed Brown 167c4762a1bSJed Brown test: 168c4762a1bSJed Brown suffix: 4_baij 169c4762a1bSJed Brown nsize: 3 170c4762a1bSJed Brown args: -mat_type baij -nsubcomms 2 -view_mat -psubcomm_type interlaced 171c4762a1bSJed Brown 172c4762a1bSJed Brown test: 173c4762a1bSJed Brown suffix: 4_sbaij 174c4762a1bSJed Brown nsize: 3 175c4762a1bSJed Brown args: -mat_type sbaij -nsubcomms 2 -view_mat -psubcomm_type interlaced 176c4762a1bSJed Brown 177c4762a1bSJed Brown test: 178c4762a1bSJed Brown suffix: 4_dense 179c4762a1bSJed Brown nsize: 3 180c4762a1bSJed Brown args: -mat_type dense -nsubcomms 2 -view_mat -psubcomm_type interlaced 181c4762a1bSJed Brown 182c4762a1bSJed Brown TEST*/ 183