static char help[] = "Tests MatConvert(), MatLoad() for MATSCALAPACK interface.\n\n"; /* Example: mpiexec -n ./ex244 -fA -fB -orig_mat_type -orig_mat_type */ #include int main(int argc,char **args) { Mat A,Ae,B,Be; PetscViewer view; char file[2][PETSC_MAX_PATH_LEN]; PetscBool flg,flgB,isScaLAPACK,isDense,isAij,isSbaij; PetscScalar one = 1.0; PetscMPIInt rank,size; PetscInt M,N; PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); /* Load PETSc matrices */ PetscCall(PetscOptionsGetString(NULL,NULL,"-fA",file[0],PETSC_MAX_PATH_LEN,NULL)); PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[0],FILE_MODE_READ,&view)); PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); PetscCall(MatSetOptionsPrefix(A,"orig_")); PetscCall(MatSetType(A,MATAIJ)); PetscCall(MatSetFromOptions(A)); PetscCall(MatLoad(A,view)); PetscCall(PetscViewerDestroy(&view)); PetscOptionsGetString(NULL,NULL,"-fB",file[1],PETSC_MAX_PATH_LEN,&flgB); if (flgB) { PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file[1],FILE_MODE_READ,&view)); PetscCall(MatCreate(PETSC_COMM_WORLD,&B)); PetscCall(MatSetOptionsPrefix(B,"orig_")); PetscCall(MatSetType(B,MATAIJ)); PetscCall(MatSetFromOptions(B)); PetscCall(MatLoad(B,view)); PetscCall(PetscViewerDestroy(&view)); } else { /* Create matrix B = I */ PetscInt rstart,rend,i; PetscCall(MatGetSize(A,&M,&N)); PetscCall(MatGetOwnershipRange(A,&rstart,&rend)); PetscCall(MatCreate(PETSC_COMM_WORLD,&B)); PetscCall(MatSetOptionsPrefix(B,"orig_")); PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,M,N)); PetscCall(MatSetType(B,MATAIJ)); PetscCall(MatSetFromOptions(B)); PetscCall(MatSetUp(B)); for (i=rstart; i