1d24d4204SJose E. Roman static char help[] = "Tests ScaLAPACK interface.\n\n"; 2d24d4204SJose E. Roman 3d24d4204SJose E. Roman #include <petscmat.h> 4d24d4204SJose E. Roman 5d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 6d71ae5a4SJacob Faibussowitsch { 7f7ec113fSDamian Marek Mat Cdense, Caij, B, C, Ct, Asub; 8d24d4204SJose E. Roman Vec d; 9d24d4204SJose E. Roman PetscInt i, j, M = 5, N, mb = 2, nb, nrows, ncols, mloc, nloc; 10d24d4204SJose E. Roman const PetscInt *rows, *cols; 11d24d4204SJose E. Roman IS isrows, iscols; 12d24d4204SJose E. Roman PetscScalar *v; 13f7ec113fSDamian Marek PetscMPIInt rank, color; 14d24d4204SJose E. Roman PetscReal Cnorm; 15d24d4204SJose E. Roman PetscBool flg, mats_view = PETSC_FALSE; 16f7ec113fSDamian Marek MPI_Comm subcomm; 17d24d4204SJose E. Roman 18327415f7SBarry Smith PetscFunctionBeginUser; 19c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help)); 209566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 219566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-M", &M, NULL)); 22d24d4204SJose E. Roman N = M; 239566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-N", &N, NULL)); 249566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-mb", &mb, NULL)); 25d24d4204SJose E. Roman nb = mb; 269566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-nb", &nb, NULL)); 279566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-mats_view", &mats_view)); 28d24d4204SJose E. Roman 299566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &C)); 309566063dSJacob Faibussowitsch PetscCall(MatSetType(C, MATSCALAPACK)); 31d24d4204SJose E. Roman mloc = PETSC_DECIDE; 329566063dSJacob Faibussowitsch PetscCall(PetscSplitOwnershipEqual(PETSC_COMM_WORLD, &mloc, &M)); 33d24d4204SJose E. Roman nloc = PETSC_DECIDE; 349566063dSJacob Faibussowitsch PetscCall(PetscSplitOwnershipEqual(PETSC_COMM_WORLD, &nloc, &N)); 359566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, mloc, nloc, M, N)); 369566063dSJacob Faibussowitsch PetscCall(MatScaLAPACKSetBlockSizes(C, mb, nb)); 379566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C)); 389566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 399566063dSJacob Faibussowitsch /*PetscCall(MatCreateScaLAPACK(PETSC_COMM_WORLD,mb,nb,M,N,0,0,&C)); */ 40d24d4204SJose E. Roman 419566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipIS(C, &isrows, &iscols)); 429566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isrows, &nrows)); 439566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isrows, &rows)); 449566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(iscols, &ncols)); 459566063dSJacob Faibussowitsch PetscCall(ISGetIndices(iscols, &cols)); 469566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrows * ncols, &v)); 47d24d4204SJose E. Roman for (i = 0; i < nrows; i++) { 48d24d4204SJose E. Roman for (j = 0; j < ncols; j++) v[i * ncols + j] = (PetscReal)(rows[i] + 1 + (cols[j] + 1) * 0.01); 49d24d4204SJose E. Roman } 509566063dSJacob Faibussowitsch PetscCall(MatSetValues(C, nrows, rows, ncols, cols, v, INSERT_VALUES)); 519566063dSJacob Faibussowitsch PetscCall(PetscFree(v)); 529566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isrows, &rows)); 539566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(iscols, &cols)); 549566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 559566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 569566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrows)); 579566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscols)); 58d24d4204SJose E. Roman 59d24d4204SJose E. Roman /* Test MatView(), MatDuplicate() and out-of-place MatConvert() */ 609566063dSJacob Faibussowitsch PetscCall(MatDuplicate(C, MAT_COPY_VALUES, &B)); 61d24d4204SJose E. Roman if (mats_view) { 629566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Duplicated C:\n")); 639566063dSJacob Faibussowitsch PetscCall(MatView(B, PETSC_VIEWER_STDOUT_WORLD)); 64d24d4204SJose E. Roman } 659566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 669566063dSJacob Faibussowitsch PetscCall(MatConvert(C, MATDENSE, MAT_INITIAL_MATRIX, &Cdense)); 679566063dSJacob Faibussowitsch PetscCall(MatMultEqual(C, Cdense, 5, &flg)); 6828b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Check fails: Cdense != C"); 69d24d4204SJose E. Roman 70d24d4204SJose E. Roman /* Test MatNorm() */ 719566063dSJacob Faibussowitsch PetscCall(MatNorm(C, NORM_1, &Cnorm)); 72d24d4204SJose E. Roman 73d24d4204SJose E. Roman /* Test MatTranspose(), MatZeroEntries() and MatGetDiagonal() */ 749566063dSJacob Faibussowitsch PetscCall(MatTranspose(C, MAT_INITIAL_MATRIX, &Ct)); 759566063dSJacob Faibussowitsch PetscCall(MatConjugate(Ct)); 76d24d4204SJose E. Roman if (mats_view) { 779566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "C's Transpose Conjugate:\n")); 789566063dSJacob Faibussowitsch PetscCall(MatView(Ct, PETSC_VIEWER_STDOUT_WORLD)); 79d24d4204SJose E. Roman } 809566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(Ct)); 819566063dSJacob Faibussowitsch if (M > N) PetscCall(MatCreateVecs(C, &d, NULL)); 829566063dSJacob Faibussowitsch else PetscCall(MatCreateVecs(C, NULL, &d)); 839566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(C, d)); 84d24d4204SJose E. Roman if (mats_view) { 859566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Diagonal of C:\n")); 869566063dSJacob Faibussowitsch PetscCall(VecView(d, PETSC_VIEWER_STDOUT_WORLD)); 87d24d4204SJose E. Roman } 88d24d4204SJose E. Roman if (M > N) { 899566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(C, NULL, d)); 90d24d4204SJose E. Roman } else { 919566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(C, d, NULL)); 92d24d4204SJose E. Roman } 93d24d4204SJose E. Roman if (mats_view) { 949566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Diagonal Scaled C:\n")); 959566063dSJacob Faibussowitsch PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD)); 96d24d4204SJose E. Roman } 97d24d4204SJose E. Roman 98d24d4204SJose E. Roman /* Test MatAXPY(), MatAYPX() and in-place MatConvert() */ 999566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &B)); 1009566063dSJacob Faibussowitsch PetscCall(MatSetType(B, MATSCALAPACK)); 1019566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B, mloc, nloc, PETSC_DECIDE, PETSC_DECIDE)); 1029566063dSJacob Faibussowitsch PetscCall(MatScaLAPACKSetBlockSizes(B, mb, nb)); 1039566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(B)); 1049566063dSJacob Faibussowitsch PetscCall(MatSetUp(B)); 1059566063dSJacob Faibussowitsch /* PetscCall(MatCreateScaLAPACK(PETSC_COMM_WORLD,mb,nb,M,N,0,0,&B)); */ 1069566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipIS(B, &isrows, &iscols)); 1079566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isrows, &nrows)); 1089566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isrows, &rows)); 1099566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(iscols, &ncols)); 1109566063dSJacob Faibussowitsch PetscCall(ISGetIndices(iscols, &cols)); 1119566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrows * ncols, &v)); 112d24d4204SJose E. Roman for (i = 0; i < nrows; i++) { 113d24d4204SJose E. Roman for (j = 0; j < ncols; j++) v[i * ncols + j] = (PetscReal)(1000 * rows[i] + cols[j]); 114d24d4204SJose E. Roman } 1159566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, nrows, rows, ncols, cols, v, INSERT_VALUES)); 1169566063dSJacob Faibussowitsch PetscCall(PetscFree(v)); 1179566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isrows, &rows)); 1189566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(iscols, &cols)); 1199566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 1209566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 121d24d4204SJose E. Roman if (mats_view) { 1229566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "B:\n")); 1239566063dSJacob Faibussowitsch PetscCall(MatView(B, PETSC_VIEWER_STDOUT_WORLD)); 124d24d4204SJose E. Roman } 1259566063dSJacob Faibussowitsch PetscCall(MatAXPY(B, 2.5, C, SAME_NONZERO_PATTERN)); 1269566063dSJacob Faibussowitsch PetscCall(MatAYPX(B, 3.75, C, SAME_NONZERO_PATTERN)); 1279566063dSJacob Faibussowitsch PetscCall(MatConvert(B, MATDENSE, MAT_INPLACE_MATRIX, &B)); 128d24d4204SJose E. Roman if (mats_view) { 1299566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "B after MatAXPY and MatAYPX:\n")); 1309566063dSJacob Faibussowitsch PetscCall(MatView(B, PETSC_VIEWER_STDOUT_WORLD)); 131d24d4204SJose E. Roman } 1329566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrows)); 1339566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscols)); 1349566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 135d24d4204SJose E. Roman 136d24d4204SJose E. Roman /* Test MatMatTransposeMult(): B = C*C^T */ 137fb842aefSJose E. Roman PetscCall(MatMatTransposeMult(C, C, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &B)); 1389566063dSJacob Faibussowitsch PetscCall(MatScale(C, 2.0)); 139fb842aefSJose E. Roman PetscCall(MatMatTransposeMult(C, C, MAT_REUSE_MATRIX, PETSC_DETERMINE, &B)); 1409566063dSJacob Faibussowitsch PetscCall(MatMatTransposeMultEqual(C, C, B, 10, &flg)); 14128b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Check fails: B != C*C^T"); 142d24d4204SJose E. Roman 143d24d4204SJose E. Roman if (mats_view) { 1449566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "C MatMatTransposeMult C:\n")); 1459566063dSJacob Faibussowitsch PetscCall(MatView(B, PETSC_VIEWER_STDOUT_WORLD)); 146d24d4204SJose E. Roman } 147e0d15688SBlanca Mellado Pinto PetscCall(MatDestroy(&B)); 148e0d15688SBlanca Mellado Pinto 149e0d15688SBlanca Mellado Pinto /* Test MatTransposeMatMult(): B = C^T*C */ 150fb842aefSJose E. Roman PetscCall(MatTransposeMatMult(C, C, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &B)); 151e0d15688SBlanca Mellado Pinto PetscCall(MatScale(C, 2.0)); 152fb842aefSJose E. Roman PetscCall(MatTransposeMatMult(C, C, MAT_REUSE_MATRIX, PETSC_DETERMINE, &B)); 153e0d15688SBlanca Mellado Pinto PetscCall(MatTransposeMatMultEqual(C, C, B, 10, &flg)); 154e0d15688SBlanca Mellado Pinto PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Check fails: B != C^T*C"); 155e0d15688SBlanca Mellado Pinto 156e0d15688SBlanca Mellado Pinto if (mats_view) { 157e0d15688SBlanca Mellado Pinto PetscCall(PetscPrintf(PETSC_COMM_WORLD, "C MatTransposeMatMult C:\n")); 158e0d15688SBlanca Mellado Pinto PetscCall(MatView(B, PETSC_VIEWER_STDOUT_WORLD)); 159e0d15688SBlanca Mellado Pinto } 160d24d4204SJose E. Roman 161d24d4204SJose E. Roman /* Test MatMult() */ 1629566063dSJacob Faibussowitsch PetscCall(MatComputeOperator(C, MATAIJ, &Caij)); 1639566063dSJacob Faibussowitsch PetscCall(MatMultEqual(C, Caij, 5, &flg)); 16428b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_ARG_NOTSAMETYPE, "C != Caij. MatMultEqual() fails"); 1659566063dSJacob Faibussowitsch PetscCall(MatMultTransposeEqual(C, Caij, 5, &flg)); 16628b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_ARG_NOTSAMETYPE, "C != Caij. MatMultTransposeEqual() fails"); 167d24d4204SJose E. Roman 168d24d4204SJose E. Roman /* Test MatMultAdd() and MatMultTransposeAddEqual() */ 1699566063dSJacob Faibussowitsch PetscCall(MatMultAddEqual(C, Caij, 5, &flg)); 17028b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_ARG_NOTSAMETYPE, "C != Caij. MatMultAddEqual() fails"); 1719566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAddEqual(C, Caij, 5, &flg)); 17228b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_ARG_NOTSAMETYPE, "C != Caij. MatMultTransposeAddEqual() fails"); 173d24d4204SJose E. Roman 174d24d4204SJose E. Roman /* Test MatMatMult() */ 1759566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-test_matmatmult", &flg)); 176d24d4204SJose E. Roman if (flg) { 177d24d4204SJose E. Roman Mat CC, CCaij; 178fb842aefSJose E. Roman PetscCall(MatMatMult(C, C, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &CC)); 179fb842aefSJose E. Roman PetscCall(MatMatMult(Caij, Caij, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &CCaij)); 1809566063dSJacob Faibussowitsch PetscCall(MatMultEqual(CC, CCaij, 5, &flg)); 18128b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_ARG_NOTSAMETYPE, "CC != CCaij. MatMatMult() fails"); 1829566063dSJacob Faibussowitsch PetscCall(MatDestroy(&CCaij)); 1839566063dSJacob Faibussowitsch PetscCall(MatDestroy(&CC)); 184d24d4204SJose E. Roman } 185d24d4204SJose E. Roman 186f7ec113fSDamian Marek /* Test MatCreate() on subcomm */ 187f7ec113fSDamian Marek color = rank % 2; 1889566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_split(PETSC_COMM_WORLD, color, 0, &subcomm)); 189f7ec113fSDamian Marek if (color == 0) { 1909566063dSJacob Faibussowitsch PetscCall(MatCreate(subcomm, &Asub)); 1919566063dSJacob Faibussowitsch PetscCall(MatSetType(Asub, MATSCALAPACK)); 192f7ec113fSDamian Marek mloc = PETSC_DECIDE; 1939566063dSJacob Faibussowitsch PetscCall(PetscSplitOwnershipEqual(subcomm, &mloc, &M)); 194f7ec113fSDamian Marek nloc = PETSC_DECIDE; 1959566063dSJacob Faibussowitsch PetscCall(PetscSplitOwnershipEqual(subcomm, &nloc, &N)); 1969566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Asub, mloc, nloc, M, N)); 1979566063dSJacob Faibussowitsch PetscCall(MatScaLAPACKSetBlockSizes(Asub, mb, nb)); 1989566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(Asub)); 1999566063dSJacob Faibussowitsch PetscCall(MatSetUp(Asub)); 2009566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Asub)); 201f7ec113fSDamian Marek } 202f7ec113fSDamian Marek 2039566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Cdense)); 2049566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Caij)); 2059566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 2069566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 2079566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Ct)); 2089566063dSJacob Faibussowitsch PetscCall(VecDestroy(&d)); 2099566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&subcomm)); 2109566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 211b122ec5aSJacob Faibussowitsch return 0; 212d24d4204SJose E. Roman } 213d24d4204SJose E. Roman 214d24d4204SJose E. Roman /*TEST 215d24d4204SJose E. Roman 216d24d4204SJose E. Roman build: 217*cf053153SJunchao Zhang requires: scalapack double 218d24d4204SJose E. Roman 219*cf053153SJunchao Zhang testset: 2203886731fSPierre Jolivet output_file: output/empty.out 221d24d4204SJose E. Roman 222d24d4204SJose E. Roman test: 223*cf053153SJunchao Zhang nsize: 2 224*cf053153SJunchao Zhang args: -mb 5 -nb 5 -M 12 -N 10 225*cf053153SJunchao Zhang 226*cf053153SJunchao Zhang test: 227d24d4204SJose E. Roman suffix: 2 228d24d4204SJose E. Roman nsize: 6 229d24d4204SJose E. Roman args: -mb 8 -nb 6 -M 20 -N 50 230d24d4204SJose E. Roman 231d24d4204SJose E. Roman test: 232d24d4204SJose E. Roman suffix: 3 233d24d4204SJose E. Roman nsize: 3 234d24d4204SJose E. Roman args: -mb 2 -nb 2 -M 20 -N 20 -test_matmatmult 235d24d4204SJose E. Roman 236d24d4204SJose E. Roman TEST*/ 237