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