1 2 static char help[] = "Tests MatConvert(), MatLoad() for MATSCALAPACK interface.\n\n"; 3 /* 4 Example: 5 mpiexec -n <np> ./ex244 -fA <A_data> -fB <B_data> -orig_mat_type <type> -orig_mat_type <mat_type> 6 */ 7 8 #include <petscmat.h> 9 10 int main(int argc,char **args) 11 { 12 Mat A,Ae,B,Be; 13 PetscViewer view; 14 char file[2][PETSC_MAX_PATH_LEN]; 15 PetscBool flg,flgB,isScaLAPACK,isDense,isAij,isSbaij; 16 PetscScalar one = 1.0; 17 PetscMPIInt rank,size; 18 PetscInt M,N; 19 20 PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 21 PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 22 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 23 24 /* Load PETSc matrices */ 25 PetscCall(PetscOptionsGetString(NULL,NULL,"-fA",file[0],PETSC_MAX_PATH_LEN,NULL)); 26 PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[0],FILE_MODE_READ,&view)); 27 PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); 28 PetscCall(MatSetOptionsPrefix(A,"orig_")); 29 PetscCall(MatSetType(A,MATAIJ)); 30 PetscCall(MatSetFromOptions(A)); 31 PetscCall(MatLoad(A,view)); 32 PetscCall(PetscViewerDestroy(&view)); 33 34 PetscOptionsGetString(NULL,NULL,"-fB",file[1],PETSC_MAX_PATH_LEN,&flgB); 35 if (flgB) { 36 PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[1],FILE_MODE_READ,&view)); 37 PetscCall(MatCreate(PETSC_COMM_WORLD,&B)); 38 PetscCall(MatSetOptionsPrefix(B,"orig_")); 39 PetscCall(MatSetType(B,MATAIJ)); 40 PetscCall(MatSetFromOptions(B)); 41 PetscCall(MatLoad(B,view)); 42 PetscCall(PetscViewerDestroy(&view)); 43 } else { 44 /* Create matrix B = I */ 45 PetscInt rstart,rend,i; 46 PetscCall(MatGetSize(A,&M,&N)); 47 PetscCall(MatGetOwnershipRange(A,&rstart,&rend)); 48 49 PetscCall(MatCreate(PETSC_COMM_WORLD,&B)); 50 PetscCall(MatSetOptionsPrefix(B,"orig_")); 51 PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,M,N)); 52 PetscCall(MatSetType(B,MATAIJ)); 53 PetscCall(MatSetFromOptions(B)); 54 PetscCall(MatSetUp(B)); 55 for (i=rstart; i<rend; i++) { 56 PetscCall(MatSetValues(B,1,&i,1,&i,&one,ADD_VALUES)); 57 } 58 PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 59 PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 60 } 61 62 PetscCall(PetscObjectTypeCompare((PetscObject)A,MATSCALAPACK,&isScaLAPACK)); 63 if (isScaLAPACK) { 64 Ae = A; 65 Be = B; 66 isDense = isAij = isSbaij = PETSC_FALSE; 67 } else { /* Convert AIJ/DENSE/SBAIJ matrices into ScaLAPACK matrices */ 68 if (size == 1) { 69 PetscCall(PetscObjectTypeCompare((PetscObject)A,MATSEQDENSE,&isDense)); 70 PetscCall(PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isAij)); 71 PetscCall(PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSbaij)); 72 } else { 73 PetscCall(PetscObjectTypeCompare((PetscObject)A,MATMPIDENSE,&isDense)); 74 PetscCall(PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isAij)); 75 PetscCall(PetscObjectTypeCompare((PetscObject)A,MATMPISBAIJ,&isSbaij)); 76 } 77 78 if (rank == 0) { 79 if (isDense) { 80 printf(" Convert DENSE matrices A and B into ScaLAPACK matrix... \n"); 81 } else if (isAij) { 82 printf(" Convert AIJ matrices A and B into ScaLAPACK matrix... \n"); 83 } else if (isSbaij) { 84 printf(" Convert SBAIJ matrices A and B into ScaLAPACK matrix... \n"); 85 } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Not supported yet"); 86 } 87 PetscCall(MatConvert(A, MATSCALAPACK, MAT_INITIAL_MATRIX, &Ae)); 88 PetscCall(MatConvert(B, MATSCALAPACK, MAT_INITIAL_MATRIX, &Be)); 89 90 /* Test accuracy */ 91 PetscCall(MatMultEqual(A,Ae,5,&flg)); 92 PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"A != A_scalapack."); 93 PetscCall(MatMultEqual(B,Be,5,&flg)); 94 PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"B != B_scalapack."); 95 } 96 97 if (!isScaLAPACK) { 98 PetscCall(MatDestroy(&Ae)); 99 PetscCall(MatDestroy(&Be)); 100 101 /* Test MAT_REUSE_MATRIX which is only supported for inplace conversion */ 102 PetscCall(MatConvert(A, MATSCALAPACK, MAT_INPLACE_MATRIX, &A)); 103 //PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 104 } 105 106 PetscCall(MatDestroy(&A)); 107 PetscCall(MatDestroy(&B)); 108 PetscCall(PetscFinalize()); 109 return 0; 110 } 111 112 /*TEST 113 114 build: 115 requires: scalapack 116 117 test: 118 requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 119 args: -fA ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_A_aij -fB ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_B_aij 120 output_file: output/ex244.out 121 122 test: 123 suffix: 2 124 nsize: 8 125 requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 126 args: -fA ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_A_aij -fB ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_B_aij 127 output_file: output/ex244.out 128 129 test: 130 suffix: 2_dense 131 nsize: 8 132 requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 133 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 134 output_file: output/ex244_dense.out 135 136 test: 137 suffix: 2_scalapack 138 nsize: 8 139 requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 140 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 141 output_file: output/ex244_scalapack.out 142 143 test: 144 suffix: 2_sbaij 145 nsize: 8 146 requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 147 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 148 output_file: output/ex244_sbaij.out 149 150 test: 151 suffix: complex 152 requires: complex double datafilespath !defined(PETSC_USE_64BIT_INDICES) 153 args: -fA ${DATAFILESPATH}/matrices/nimrod/small_112905 154 output_file: output/ex244.out 155 156 test: 157 suffix: complex_2 158 nsize: 4 159 requires: complex double datafilespath !defined(PETSC_USE_64BIT_INDICES) 160 args: -fA ${DATAFILESPATH}/matrices/nimrod/small_112905 161 output_file: output/ex244.out 162 163 test: 164 suffix: dense 165 requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 166 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 167 168 test: 169 suffix: scalapack 170 requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 171 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 172 173 test: 174 suffix: sbaij 175 requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 176 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 177 178 TEST*/ 179