1d24d4204SJose E. Roman 2d24d4204SJose E. Roman static char help[] = "Tests MatConvert(), MatLoad() for MATSCALAPACK interface.\n\n"; 3d24d4204SJose E. Roman /* 4d24d4204SJose E. Roman Example: 5d24d4204SJose E. Roman mpiexec -n <np> ./ex244 -fA <A_data> -fB <B_data> -orig_mat_type <type> -orig_mat_type <mat_type> 6d24d4204SJose E. Roman */ 7d24d4204SJose E. Roman 8d24d4204SJose E. Roman #include <petscmat.h> 9d24d4204SJose E. Roman 10d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 11d71ae5a4SJacob Faibussowitsch { 12d24d4204SJose E. Roman Mat A, Ae, B, Be; 13d24d4204SJose E. Roman PetscViewer view; 14d24d4204SJose E. Roman char file[2][PETSC_MAX_PATH_LEN]; 15d24d4204SJose E. Roman PetscBool flg, flgB, isScaLAPACK, isDense, isAij, isSbaij; 16d24d4204SJose E. Roman PetscScalar one = 1.0; 17d24d4204SJose E. Roman PetscMPIInt rank, size; 18d24d4204SJose E. Roman PetscInt M, N; 19d24d4204SJose E. Roman 20327415f7SBarry Smith PetscFunctionBeginUser; 219566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 239566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 24d24d4204SJose E. Roman 25d24d4204SJose E. Roman /* Load PETSc matrices */ 269566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-fA", file[0], PETSC_MAX_PATH_LEN, NULL)); 279566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file[0], FILE_MODE_READ, &view)); 289566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 299566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(A, "orig_")); 309566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATAIJ)); 319566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 329566063dSJacob Faibussowitsch PetscCall(MatLoad(A, view)); 339566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&view)); 34d24d4204SJose E. Roman 35*3ba16761SJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-fB", file[1], PETSC_MAX_PATH_LEN, &flgB)); 36d24d4204SJose E. Roman if (flgB) { 379566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file[1], FILE_MODE_READ, &view)); 389566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &B)); 399566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(B, "orig_")); 409566063dSJacob Faibussowitsch PetscCall(MatSetType(B, MATAIJ)); 419566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(B)); 429566063dSJacob Faibussowitsch PetscCall(MatLoad(B, view)); 439566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&view)); 44d24d4204SJose E. Roman } else { 45d24d4204SJose E. Roman /* Create matrix B = I */ 46d24d4204SJose E. Roman PetscInt rstart, rend, i; 479566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M, &N)); 489566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &rstart, &rend)); 49d24d4204SJose E. Roman 509566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &B)); 519566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(B, "orig_")); 529566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B, PETSC_DECIDE, PETSC_DECIDE, M, N)); 539566063dSJacob Faibussowitsch PetscCall(MatSetType(B, MATAIJ)); 549566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(B)); 559566063dSJacob Faibussowitsch PetscCall(MatSetUp(B)); 5648a46eb9SPierre Jolivet for (i = rstart; i < rend; i++) PetscCall(MatSetValues(B, 1, &i, 1, &i, &one, ADD_VALUES)); 579566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 589566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 59d24d4204SJose E. Roman } 60d24d4204SJose E. Roman 619566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSCALAPACK, &isScaLAPACK)); 62d24d4204SJose E. Roman if (isScaLAPACK) { 63d24d4204SJose E. Roman Ae = A; 64d24d4204SJose E. Roman Be = B; 65d24d4204SJose E. Roman isDense = isAij = isSbaij = PETSC_FALSE; 66d24d4204SJose E. Roman } else { /* Convert AIJ/DENSE/SBAIJ matrices into ScaLAPACK matrices */ 67d24d4204SJose E. Roman if (size == 1) { 689566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQDENSE, &isDense)); 699566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQAIJ, &isAij)); 709566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQSBAIJ, &isSbaij)); 71d24d4204SJose E. Roman } else { 729566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPIDENSE, &isDense)); 739566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPIAIJ, &isAij)); 749566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPISBAIJ, &isSbaij)); 75d24d4204SJose E. Roman } 76d24d4204SJose E. Roman 77dd400576SPatrick Sanan if (rank == 0) { 78d24d4204SJose E. Roman if (isDense) { 79d24d4204SJose E. Roman printf(" Convert DENSE matrices A and B into ScaLAPACK matrix... \n"); 80d24d4204SJose E. Roman } else if (isAij) { 81d24d4204SJose E. Roman printf(" Convert AIJ matrices A and B into ScaLAPACK matrix... \n"); 82d24d4204SJose E. Roman } else if (isSbaij) { 83d24d4204SJose E. Roman printf(" Convert SBAIJ matrices A and B into ScaLAPACK matrix... \n"); 84d24d4204SJose E. Roman } else SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Not supported yet"); 85d24d4204SJose E. Roman } 869566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATSCALAPACK, MAT_INITIAL_MATRIX, &Ae)); 879566063dSJacob Faibussowitsch PetscCall(MatConvert(B, MATSCALAPACK, MAT_INITIAL_MATRIX, &Be)); 88d24d4204SJose E. Roman 89d24d4204SJose E. Roman /* Test accuracy */ 909566063dSJacob Faibussowitsch PetscCall(MatMultEqual(A, Ae, 5, &flg)); 9128b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_NOTSAMETYPE, "A != A_scalapack."); 929566063dSJacob Faibussowitsch PetscCall(MatMultEqual(B, Be, 5, &flg)); 9328b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_ARG_NOTSAMETYPE, "B != B_scalapack."); 94d24d4204SJose E. Roman } 95d24d4204SJose E. Roman 96d24d4204SJose E. Roman if (!isScaLAPACK) { 979566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Ae)); 989566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Be)); 99d24d4204SJose E. Roman 100d24d4204SJose E. Roman /* Test MAT_REUSE_MATRIX which is only supported for inplace conversion */ 1019566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATSCALAPACK, MAT_INPLACE_MATRIX, &A)); 1029566063dSJacob Faibussowitsch //PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 103d24d4204SJose E. Roman } 104d24d4204SJose E. Roman 1059566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1069566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 1079566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 108b122ec5aSJacob Faibussowitsch return 0; 109d24d4204SJose E. Roman } 110d24d4204SJose E. Roman 111d24d4204SJose E. Roman /*TEST 112d24d4204SJose E. Roman 113d24d4204SJose E. Roman build: 114d24d4204SJose E. Roman requires: scalapack 115d24d4204SJose E. Roman 116d24d4204SJose E. Roman test: 117dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 118d24d4204SJose E. Roman args: -fA ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_A_aij -fB ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_B_aij 119d24d4204SJose E. Roman output_file: output/ex244.out 120d24d4204SJose E. Roman 121d24d4204SJose E. Roman test: 122d24d4204SJose E. Roman suffix: 2 123d24d4204SJose E. Roman nsize: 8 124dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 125d24d4204SJose E. Roman args: -fA ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_A_aij -fB ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_B_aij 126d24d4204SJose E. Roman output_file: output/ex244.out 127d24d4204SJose E. Roman 128d24d4204SJose E. Roman test: 129d24d4204SJose E. Roman suffix: 2_dense 130d24d4204SJose E. Roman nsize: 8 131dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 132d24d4204SJose E. Roman args: -fA ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_A_aij -fB ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_B_aij -orig_mat_type dense 133d24d4204SJose E. Roman output_file: output/ex244_dense.out 134d24d4204SJose E. Roman 135d24d4204SJose E. Roman test: 136d24d4204SJose E. Roman suffix: 2_scalapack 137d24d4204SJose E. Roman nsize: 8 138dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 139d24d4204SJose E. Roman args: -fA ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_A_aij -fB ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_B_aij -orig_mat_type scalapack 140d24d4204SJose E. Roman output_file: output/ex244_scalapack.out 141d24d4204SJose E. Roman 142d24d4204SJose E. Roman test: 143d24d4204SJose E. Roman suffix: 2_sbaij 144d24d4204SJose E. Roman nsize: 8 145dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 146d24d4204SJose E. Roman args: -fA ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_A -fB ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_B -orig_mat_type sbaij 147d24d4204SJose E. Roman output_file: output/ex244_sbaij.out 148d24d4204SJose E. Roman 149d24d4204SJose E. Roman test: 150d24d4204SJose E. Roman suffix: complex 151dfd57a17SPierre Jolivet requires: complex double datafilespath !defined(PETSC_USE_64BIT_INDICES) 152d24d4204SJose E. Roman args: -fA ${DATAFILESPATH}/matrices/nimrod/small_112905 153d24d4204SJose E. Roman output_file: output/ex244.out 154d24d4204SJose E. Roman 155d24d4204SJose E. Roman test: 156d24d4204SJose E. Roman suffix: complex_2 157d24d4204SJose E. Roman nsize: 4 158dfd57a17SPierre Jolivet requires: complex double datafilespath !defined(PETSC_USE_64BIT_INDICES) 159d24d4204SJose E. Roman args: -fA ${DATAFILESPATH}/matrices/nimrod/small_112905 160d24d4204SJose E. Roman output_file: output/ex244.out 161d24d4204SJose E. Roman 162d24d4204SJose E. Roman test: 163d24d4204SJose E. Roman suffix: dense 164dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 165d24d4204SJose E. Roman args: -fA ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_A_aij -fB ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_B_aij -orig_mat_type dense 166d24d4204SJose E. Roman 167d24d4204SJose E. Roman test: 168d24d4204SJose E. Roman suffix: scalapack 169dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 170d24d4204SJose E. Roman args: -fA ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_A_aij -fB ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_B_aij -orig_mat_type scalapack 171d24d4204SJose E. Roman 172d24d4204SJose E. Roman test: 173d24d4204SJose E. Roman suffix: sbaij 174dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 175d24d4204SJose E. Roman args: -fA ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_A -fB ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_B -orig_mat_type sbaij 176d24d4204SJose E. Roman 177d24d4204SJose E. Roman TEST*/ 178