1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests MatConvert(), MatLoad() for MATELEMENTAL interface.\n\n"; 3c4762a1bSJed Brown /* 4c4762a1bSJed Brown Example: 5c4762a1bSJed Brown mpiexec -n <np> ./ex173 -fA <A_data> -fB <B_data> -orig_mat_type <type> -orig_mat_type <mat_type> 6c4762a1bSJed Brown */ 7c4762a1bSJed Brown 8c4762a1bSJed Brown #include <petscmat.h> 9c4762a1bSJed Brown #include <petscmatelemental.h> 10c4762a1bSJed Brown 11c4762a1bSJed Brown int main(int argc,char **args) 12c4762a1bSJed Brown { 13c4762a1bSJed Brown Mat A,Ae,B,Be; 14c4762a1bSJed Brown PetscViewer view; 15c4762a1bSJed Brown char file[2][PETSC_MAX_PATH_LEN]; 16c4762a1bSJed Brown PetscBool flg,flgB,isElemental,isDense,isAij,isSbaij; 17c4762a1bSJed Brown PetscScalar one = 1.0; 18c4762a1bSJed Brown PetscMPIInt rank,size; 19c4762a1bSJed Brown PetscInt M,N; 20c4762a1bSJed Brown 21*327415f7SBarry Smith PetscFunctionBeginUser; 229566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 239566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 249566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 25c4762a1bSJed Brown 26c4762a1bSJed Brown /* Load PETSc matrices */ 279566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL,NULL,"-fA",file[0],sizeof(file[0]),NULL)); 289566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[0],FILE_MODE_READ,&view)); 299566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); 309566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(A,"orig_")); 319566063dSJacob Faibussowitsch PetscCall(MatSetType(A,MATAIJ)); 329566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 339566063dSJacob Faibussowitsch PetscCall(MatLoad(A,view)); 349566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&view)); 35c4762a1bSJed Brown 36589a23caSBarry Smith PetscOptionsGetString(NULL,NULL,"-fB",file[1],sizeof(file[1]),&flgB); 37c4762a1bSJed Brown if (flgB) { 389566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[1],FILE_MODE_READ,&view)); 399566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&B)); 409566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(B,"orig_")); 419566063dSJacob Faibussowitsch PetscCall(MatSetType(B,MATAIJ)); 429566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(B)); 439566063dSJacob Faibussowitsch PetscCall(MatLoad(B,view)); 449566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&view)); 45c4762a1bSJed Brown } else { 46c4762a1bSJed Brown /* Create matrix B = I */ 47c4762a1bSJed Brown PetscInt rstart,rend,i; 489566063dSJacob Faibussowitsch PetscCall(MatGetSize(A,&M,&N)); 499566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&rstart,&rend)); 50c4762a1bSJed Brown 519566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&B)); 529566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(B,"orig_")); 539566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,M,N)); 549566063dSJacob Faibussowitsch PetscCall(MatSetType(B,MATAIJ)); 559566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(B)); 569566063dSJacob Faibussowitsch PetscCall(MatSetUp(B)); 57c4762a1bSJed Brown for (i=rstart; i<rend; i++) { 589566063dSJacob Faibussowitsch PetscCall(MatSetValues(B,1,&i,1,&i,&one,ADD_VALUES)); 59c4762a1bSJed Brown } 609566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 619566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 62c4762a1bSJed Brown } 63c4762a1bSJed Brown 649566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A,MATELEMENTAL,&isElemental)); 65c4762a1bSJed Brown if (isElemental) { 66c4762a1bSJed Brown Ae = A; 67c4762a1bSJed Brown Be = B; 68c4762a1bSJed Brown isDense = isAij = isSbaij = PETSC_FALSE; 69c4762a1bSJed Brown } else { /* Convert AIJ/DENSE/SBAIJ matrices into Elemental matrices */ 70c4762a1bSJed Brown if (size == 1) { 719566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A,MATSEQDENSE,&isDense)); 729566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isAij)); 739566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSbaij)); 74c4762a1bSJed Brown } else { 759566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A,MATMPIDENSE,&isDense)); 769566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isAij)); 779566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)A,MATMPISBAIJ,&isSbaij)); 78c4762a1bSJed Brown } 79c4762a1bSJed Brown 80dd400576SPatrick Sanan if (rank == 0) { 81c4762a1bSJed Brown if (isDense) { 82c4762a1bSJed Brown printf(" Convert DENSE matrices A and B into Elemental matrix... \n"); 83c4762a1bSJed Brown } else if (isAij) { 84c4762a1bSJed Brown printf(" Convert AIJ matrices A and B into Elemental matrix... \n"); 85c4762a1bSJed Brown } else if (isSbaij) { 86c4762a1bSJed Brown printf(" Convert SBAIJ matrices A and B into Elemental matrix... \n"); 87c4762a1bSJed Brown } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Not supported yet"); 88c4762a1bSJed Brown } 899566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATELEMENTAL, MAT_INITIAL_MATRIX, &Ae)); 909566063dSJacob Faibussowitsch PetscCall(MatConvert(B, MATELEMENTAL, MAT_INITIAL_MATRIX, &Be)); 91c4762a1bSJed Brown 92c4762a1bSJed Brown /* Test accuracy */ 939566063dSJacob Faibussowitsch PetscCall(MatMultEqual(A,Ae,5,&flg)); 9428b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"A != A_elemental."); 959566063dSJacob Faibussowitsch PetscCall(MatMultEqual(B,Be,5,&flg)); 9628b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"B != B_elemental."); 97c4762a1bSJed Brown } 98c4762a1bSJed Brown 99c4762a1bSJed Brown if (!isElemental) { 1009566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Ae)); 1019566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Be)); 102c4762a1bSJed Brown 103c4762a1bSJed Brown /* Test MAT_REUSE_MATRIX which is only supported for inplace conversion */ 1049566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATELEMENTAL, MAT_INPLACE_MATRIX, &A)); 1059566063dSJacob Faibussowitsch //PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 106c4762a1bSJed Brown } 107c4762a1bSJed Brown 1089566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1099566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 1109566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 111b122ec5aSJacob Faibussowitsch return 0; 112c4762a1bSJed Brown } 113c4762a1bSJed Brown 114c4762a1bSJed Brown /*TEST 115c4762a1bSJed Brown 116c4762a1bSJed Brown build: 117c4762a1bSJed Brown requires: elemental 118c4762a1bSJed Brown 119c4762a1bSJed Brown test: 120dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 121c4762a1bSJed Brown args: -fA ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_A_aij -fB ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_B_aij 122c4762a1bSJed Brown output_file: output/ex174.out 123c4762a1bSJed Brown 124c4762a1bSJed Brown test: 125c4762a1bSJed Brown suffix: 2 126c4762a1bSJed Brown nsize: 8 127dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 128c4762a1bSJed Brown args: -fA ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_A_aij -fB ${DATAFILESPATH}/matrices/EigenProblems/Eigdftb/dftb_bin/graphene_xxs_B_aij 129c4762a1bSJed Brown output_file: output/ex174.out 130c4762a1bSJed Brown 131c4762a1bSJed Brown test: 132c4762a1bSJed Brown suffix: 2_dense 133c4762a1bSJed Brown nsize: 8 134dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 135c4762a1bSJed Brown 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 136c4762a1bSJed Brown output_file: output/ex174_dense.out 137c4762a1bSJed Brown 138c4762a1bSJed Brown test: 139c4762a1bSJed Brown suffix: 2_elemental 140c4762a1bSJed Brown nsize: 8 141dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 142c4762a1bSJed Brown 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 elemental 143c4762a1bSJed Brown output_file: output/ex174_elemental.out 144c4762a1bSJed Brown 145c4762a1bSJed Brown test: 146c4762a1bSJed Brown suffix: 2_sbaij 147c4762a1bSJed Brown nsize: 8 148dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 149c4762a1bSJed Brown 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 150c4762a1bSJed Brown output_file: output/ex174_sbaij.out 151c4762a1bSJed Brown 152c4762a1bSJed Brown test: 153c4762a1bSJed Brown suffix: complex 154dfd57a17SPierre Jolivet requires: complex double datafilespath !defined(PETSC_USE_64BIT_INDICES) 155c4762a1bSJed Brown args: -fA ${DATAFILESPATH}/matrices/nimrod/small_112905 156c4762a1bSJed Brown output_file: output/ex174.out 157c4762a1bSJed Brown 158c4762a1bSJed Brown test: 159c4762a1bSJed Brown suffix: complex_2 160c4762a1bSJed Brown nsize: 4 161dfd57a17SPierre Jolivet requires: complex double datafilespath !defined(PETSC_USE_64BIT_INDICES) 162c4762a1bSJed Brown args: -fA ${DATAFILESPATH}/matrices/nimrod/small_112905 163c4762a1bSJed Brown output_file: output/ex174.out 164c4762a1bSJed Brown 165c4762a1bSJed Brown test: 166c4762a1bSJed Brown suffix: dense 167dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 168c4762a1bSJed Brown 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 169c4762a1bSJed Brown 170c4762a1bSJed Brown test: 171c4762a1bSJed Brown suffix: elemental 172dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 173c4762a1bSJed Brown 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 elemental 174c4762a1bSJed Brown 175c4762a1bSJed Brown test: 176c4762a1bSJed Brown suffix: sbaij 177dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 178c4762a1bSJed Brown 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 179c4762a1bSJed Brown 180c4762a1bSJed Brown TEST*/ 181