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