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