1d24d4204SJose E. Roman static char help[] = "Test conversion of ScaLAPACK matrices.\n\n"; 2d24d4204SJose E. Roman 3d24d4204SJose E. Roman #include <petscmat.h> 4d24d4204SJose E. Roman 5d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 6d71ae5a4SJacob Faibussowitsch { 7d24d4204SJose E. Roman Mat A, A_scalapack; 8d24d4204SJose E. Roman PetscInt i, j, M = 10, N = 5, nloc, mloc, nrows, ncols; 9d24d4204SJose E. Roman PetscMPIInt rank, size; 10d24d4204SJose E. Roman IS isrows, iscols; 11d24d4204SJose E. Roman const PetscInt *rows, *cols; 12d24d4204SJose E. Roman PetscScalar *v; 13d24d4204SJose E. Roman MatType type; 14d24d4204SJose E. Roman PetscBool isDense, isAIJ, flg; 15d24d4204SJose E. Roman 16327415f7SBarry Smith PetscFunctionBeginUser; 17c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 189566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 199566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 209566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-M", &M, NULL)); 219566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-N", &N, NULL)); 22d24d4204SJose E. Roman 23d24d4204SJose E. Roman /* Create a matrix */ 249566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 25d24d4204SJose E. Roman mloc = PETSC_DECIDE; 269566063dSJacob Faibussowitsch PetscCall(PetscSplitOwnershipEqual(PETSC_COMM_WORLD, &mloc, &M)); 27d24d4204SJose E. Roman nloc = PETSC_DECIDE; 289566063dSJacob Faibussowitsch PetscCall(PetscSplitOwnershipEqual(PETSC_COMM_WORLD, &nloc, &N)); 299566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, mloc, nloc, M, N)); 309566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATDENSE)); 319566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 329566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 33d24d4204SJose E. Roman 34d24d4204SJose E. Roman /* Set local matrix entries */ 359566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipIS(A, &isrows, &iscols)); 369566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isrows, &nrows)); 379566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isrows, &rows)); 389566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(iscols, &ncols)); 399566063dSJacob Faibussowitsch PetscCall(ISGetIndices(iscols, &cols)); 409566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrows * ncols, &v)); 41d24d4204SJose E. Roman 42d24d4204SJose E. Roman for (i = 0; i < nrows; i++) { 43d24d4204SJose E. Roman for (j = 0; j < ncols; j++) { 44d24d4204SJose E. Roman if (size == 1) { 45d24d4204SJose E. Roman v[i * ncols + j] = (PetscScalar)(i + j); 46d24d4204SJose E. Roman } else { 47d24d4204SJose E. Roman v[i * ncols + j] = (PetscScalar)rank + j * 0.1; 48d24d4204SJose E. Roman } 49d24d4204SJose E. Roman } 50d24d4204SJose E. Roman } 519566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, nrows, rows, ncols, cols, v, INSERT_VALUES)); 529566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 539566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 54d24d4204SJose E. Roman 55d24d4204SJose E. Roman /* Test MatSetValues() by converting A to A_scalapack */ 569566063dSJacob Faibussowitsch PetscCall(MatGetType(A, &type)); 57d24d4204SJose E. Roman if (size == 1) { 589566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQDENSE, &isDense)); 599566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQAIJ, &isAIJ)); 60d24d4204SJose E. Roman } else { 619566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPIDENSE, &isDense)); 629566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPIAIJ, &isAIJ)); 63d24d4204SJose E. Roman } 64d24d4204SJose E. Roman 65d24d4204SJose E. Roman if (isDense || isAIJ) { 66d24d4204SJose E. Roman Mat Aexplicit; 679566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATSCALAPACK, MAT_INITIAL_MATRIX, &A_scalapack)); 689566063dSJacob Faibussowitsch PetscCall(MatComputeOperator(A_scalapack, isAIJ ? MATAIJ : MATDENSE, &Aexplicit)); 699566063dSJacob Faibussowitsch PetscCall(MatMultEqual(Aexplicit, A_scalapack, 5, &flg)); 7028b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Aexplicit != A_scalapack."); 719566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Aexplicit)); 72d24d4204SJose E. Roman 73d24d4204SJose E. Roman /* Test MAT_REUSE_MATRIX which is only supported for inplace conversion */ 749566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATSCALAPACK, MAT_INPLACE_MATRIX, &A)); 759566063dSJacob Faibussowitsch PetscCall(MatMultEqual(A_scalapack, A, 5, &flg)); 7628b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "A_scalapack != A."); 779566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A_scalapack)); 78d24d4204SJose E. Roman } 79d24d4204SJose E. Roman 809566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isrows, &rows)); 819566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(iscols, &cols)); 829566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrows)); 839566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscols)); 849566063dSJacob Faibussowitsch PetscCall(PetscFree(v)); 859566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 869566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 87b122ec5aSJacob Faibussowitsch return 0; 88d24d4204SJose E. Roman } 89d24d4204SJose E. Roman 90d24d4204SJose E. Roman /*TEST 91d24d4204SJose E. Roman 92d24d4204SJose E. Roman build: 93*cf053153SJunchao Zhang requires: scalapack double 94d24d4204SJose E. Roman 95*cf053153SJunchao Zhang testset: 963886731fSPierre Jolivet output_file: output/empty.out 97*cf053153SJunchao Zhang nsize: 6 98d24d4204SJose E. Roman 99d24d4204SJose E. Roman test: 100*cf053153SJunchao Zhang 101*cf053153SJunchao Zhang test: 102d24d4204SJose E. Roman suffix: 2 103d24d4204SJose E. Roman args: -mat_type aij 104d24d4204SJose E. Roman 105d24d4204SJose E. Roman test: 106d24d4204SJose E. Roman suffix: 3 107d24d4204SJose E. Roman args: -mat_type scalapack 108d24d4204SJose E. Roman TEST*/ 109