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 10d24d4204SJose E. Roman int main(int argc,char **args) 11d24d4204SJose E. Roman { 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 20*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 215f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 225f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 23d24d4204SJose E. Roman 24d24d4204SJose E. Roman /* Load PETSc matrices */ 255f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetString(NULL,NULL,"-fA",file[0],PETSC_MAX_PATH_LEN,NULL)); 265f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[0],FILE_MODE_READ,&view)); 275f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 285f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOptionsPrefix(A,"orig_")); 295f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A,MATAIJ)); 305f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(A)); 315f80ce2aSJacob Faibussowitsch CHKERRQ(MatLoad(A,view)); 325f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&view)); 33d24d4204SJose E. Roman 34d24d4204SJose E. Roman PetscOptionsGetString(NULL,NULL,"-fB",file[1],PETSC_MAX_PATH_LEN,&flgB); 35d24d4204SJose E. Roman if (flgB) { 365f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[1],FILE_MODE_READ,&view)); 375f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&B)); 385f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOptionsPrefix(B,"orig_")); 395f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(B,MATAIJ)); 405f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(B)); 415f80ce2aSJacob Faibussowitsch CHKERRQ(MatLoad(B,view)); 425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&view)); 43d24d4204SJose E. Roman } else { 44d24d4204SJose E. Roman /* Create matrix B = I */ 45d24d4204SJose E. Roman PetscInt rstart,rend,i; 465f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(A,&M,&N)); 475f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(A,&rstart,&rend)); 48d24d4204SJose E. Roman 495f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&B)); 505f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOptionsPrefix(B,"orig_")); 515f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,M,N)); 525f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(B,MATAIJ)); 535f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(B)); 545f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(B)); 55d24d4204SJose E. Roman for (i=rstart; i<rend; i++) { 565f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(B,1,&i,1,&i,&one,ADD_VALUES)); 57d24d4204SJose E. Roman } 585f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 595f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 60d24d4204SJose E. Roman } 61d24d4204SJose E. Roman 625f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)A,MATSCALAPACK,&isScaLAPACK)); 63d24d4204SJose E. Roman if (isScaLAPACK) { 64d24d4204SJose E. Roman Ae = A; 65d24d4204SJose E. Roman Be = B; 66d24d4204SJose E. Roman isDense = isAij = isSbaij = PETSC_FALSE; 67d24d4204SJose E. Roman } else { /* Convert AIJ/DENSE/SBAIJ matrices into ScaLAPACK matrices */ 68d24d4204SJose E. Roman if (size == 1) { 695f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)A,MATSEQDENSE,&isDense)); 705f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isAij)); 715f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSbaij)); 72d24d4204SJose E. Roman } else { 735f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)A,MATMPIDENSE,&isDense)); 745f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isAij)); 755f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)A,MATMPISBAIJ,&isSbaij)); 76d24d4204SJose E. Roman } 77d24d4204SJose E. Roman 78dd400576SPatrick Sanan if (rank == 0) { 79d24d4204SJose E. Roman if (isDense) { 80d24d4204SJose E. Roman printf(" Convert DENSE matrices A and B into ScaLAPACK matrix... \n"); 81d24d4204SJose E. Roman } else if (isAij) { 82d24d4204SJose E. Roman printf(" Convert AIJ matrices A and B into ScaLAPACK matrix... \n"); 83d24d4204SJose E. Roman } else if (isSbaij) { 84d24d4204SJose E. Roman printf(" Convert SBAIJ matrices A and B into ScaLAPACK matrix... \n"); 85d24d4204SJose E. Roman } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Not supported yet"); 86d24d4204SJose E. Roman } 875f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(A, MATSCALAPACK, MAT_INITIAL_MATRIX, &Ae)); 885f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(B, MATSCALAPACK, MAT_INITIAL_MATRIX, &Be)); 89d24d4204SJose E. Roman 90d24d4204SJose E. Roman /* Test accuracy */ 915f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(A,Ae,5,&flg)); 9228b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"A != A_scalapack."); 935f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(B,Be,5,&flg)); 9428b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"B != B_scalapack."); 95d24d4204SJose E. Roman } 96d24d4204SJose E. Roman 97d24d4204SJose E. Roman if (!isScaLAPACK) { 985f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Ae)); 995f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Be)); 100d24d4204SJose E. Roman 101d24d4204SJose E. Roman /* Test MAT_REUSE_MATRIX which is only supported for inplace conversion */ 1025f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(A, MATSCALAPACK, MAT_INPLACE_MATRIX, &A)); 1035f80ce2aSJacob Faibussowitsch //CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 104d24d4204SJose E. Roman } 105d24d4204SJose E. Roman 1065f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 1075f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 108*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 109*b122ec5aSJacob Faibussowitsch return 0; 110d24d4204SJose E. Roman } 111d24d4204SJose E. Roman 112d24d4204SJose E. Roman /*TEST 113d24d4204SJose E. Roman 114d24d4204SJose E. Roman build: 115d24d4204SJose E. Roman requires: scalapack 116d24d4204SJose E. Roman 117d24d4204SJose E. Roman test: 118dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 119d24d4204SJose 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 120d24d4204SJose E. Roman output_file: output/ex244.out 121d24d4204SJose E. Roman 122d24d4204SJose E. Roman test: 123d24d4204SJose E. Roman suffix: 2 124d24d4204SJose E. Roman nsize: 8 125dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 126d24d4204SJose 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 127d24d4204SJose E. Roman output_file: output/ex244.out 128d24d4204SJose E. Roman 129d24d4204SJose E. Roman test: 130d24d4204SJose E. Roman suffix: 2_dense 131d24d4204SJose E. Roman nsize: 8 132dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 133d24d4204SJose 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 134d24d4204SJose E. Roman output_file: output/ex244_dense.out 135d24d4204SJose E. Roman 136d24d4204SJose E. Roman test: 137d24d4204SJose E. Roman suffix: 2_scalapack 138d24d4204SJose E. Roman nsize: 8 139dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 140d24d4204SJose 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 141d24d4204SJose E. Roman output_file: output/ex244_scalapack.out 142d24d4204SJose E. Roman 143d24d4204SJose E. Roman test: 144d24d4204SJose E. Roman suffix: 2_sbaij 145d24d4204SJose E. Roman nsize: 8 146dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 147d24d4204SJose 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 148d24d4204SJose E. Roman output_file: output/ex244_sbaij.out 149d24d4204SJose E. Roman 150d24d4204SJose E. Roman test: 151d24d4204SJose E. Roman suffix: complex 152dfd57a17SPierre Jolivet requires: complex double datafilespath !defined(PETSC_USE_64BIT_INDICES) 153d24d4204SJose E. Roman args: -fA ${DATAFILESPATH}/matrices/nimrod/small_112905 154d24d4204SJose E. Roman output_file: output/ex244.out 155d24d4204SJose E. Roman 156d24d4204SJose E. Roman test: 157d24d4204SJose E. Roman suffix: complex_2 158d24d4204SJose E. Roman nsize: 4 159dfd57a17SPierre Jolivet requires: complex double datafilespath !defined(PETSC_USE_64BIT_INDICES) 160d24d4204SJose E. Roman args: -fA ${DATAFILESPATH}/matrices/nimrod/small_112905 161d24d4204SJose E. Roman output_file: output/ex244.out 162d24d4204SJose E. Roman 163d24d4204SJose E. Roman test: 164d24d4204SJose E. Roman suffix: dense 165dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 166d24d4204SJose 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 167d24d4204SJose E. Roman 168d24d4204SJose E. Roman test: 169d24d4204SJose E. Roman suffix: scalapack 170dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 171d24d4204SJose 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 172d24d4204SJose E. Roman 173d24d4204SJose E. Roman test: 174d24d4204SJose E. Roman suffix: sbaij 175dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 176d24d4204SJose 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 177d24d4204SJose E. Roman 178d24d4204SJose E. Roman TEST*/ 179