1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Reads a matrix from PETSc binary file. Use for view or investigating matrix data structure. \n\n"; 3c4762a1bSJed Brown /* 4c4762a1bSJed Brown Example: 5c4762a1bSJed Brown ./ex16 -f <matrix file> -a_mat_view draw -draw_pause -1 6c4762a1bSJed Brown ./ex16 -f <matrix file> -a_mat_view ascii::ascii_info 7c4762a1bSJed Brown */ 8c4762a1bSJed Brown 9c4762a1bSJed Brown #include <petscmat.h> 10c4762a1bSJed Brown int main(int argc,char **args) 11c4762a1bSJed Brown { 12c4762a1bSJed Brown Mat A,Asp; 13c4762a1bSJed Brown PetscViewer fd; /* viewer */ 14c4762a1bSJed Brown char file[PETSC_MAX_PATH_LEN]; /* input file name */ 15c4762a1bSJed Brown PetscInt m,n,rstart,rend; 16c4762a1bSJed Brown PetscBool flg; 17c4762a1bSJed Brown PetscInt row,ncols,j,nrows,nnzA=0,nnzAsp=0; 18c4762a1bSJed Brown const PetscInt *cols; 19c4762a1bSJed Brown const PetscScalar *vals; 20c4762a1bSJed Brown PetscReal norm,percent,val,dtol=1.e-16; 21c4762a1bSJed Brown PetscMPIInt rank; 22c4762a1bSJed Brown MatInfo matinfo; 23c4762a1bSJed Brown PetscInt Dnnz,Onnz; 24c4762a1bSJed Brown 25*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 265f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 27c4762a1bSJed Brown 28c4762a1bSJed Brown /* Determine files from which we read the linear systems. */ 295f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&flg)); 3028b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f option"); 31c4762a1bSJed Brown 32c4762a1bSJed Brown /* Open binary file. Note that we use FILE_MODE_READ to indicate 33c4762a1bSJed Brown reading from this file. */ 345f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd)); 35c4762a1bSJed Brown 36c4762a1bSJed Brown /* Load the matrix; then destroy the viewer. */ 375f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 385f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOptionsPrefix(A,"a_")); 395f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(A)); 405f80ce2aSJacob Faibussowitsch CHKERRQ(MatLoad(A,fd)); 415f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&fd)); 425f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(A,&m,&n)); 435f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetInfo(A,MAT_LOCAL,&matinfo)); 44c4762a1bSJed Brown /*printf("matinfo.nz_used %g\n",matinfo.nz_used);*/ 45c4762a1bSJed Brown 46c4762a1bSJed Brown /* Get a sparse matrix Asp by dumping zero entries of A */ 475f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&Asp)); 485f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(Asp,m,n,PETSC_DECIDE,PETSC_DECIDE)); 495f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOptionsPrefix(Asp,"asp_")); 505f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(Asp)); 51c4762a1bSJed Brown Dnnz = (PetscInt)matinfo.nz_used/m + 1; 52c4762a1bSJed Brown Onnz = Dnnz/2; 53c4762a1bSJed Brown printf("Dnnz %d %d\n",Dnnz,Onnz); 545f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJSetPreallocation(Asp,Dnnz,NULL)); 555f80ce2aSJacob Faibussowitsch CHKERRQ(MatMPIAIJSetPreallocation(Asp,Dnnz,NULL,Onnz,NULL)); 56c4762a1bSJed Brown /* The allocation above is approximate so we must set this option to be permissive. 57c4762a1bSJed Brown * Real code should preallocate exactly. */ 585f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(Asp,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE)); 59c4762a1bSJed Brown 60c4762a1bSJed Brown /* Check zero rows */ 615f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(A,&rstart,&rend)); 62c4762a1bSJed Brown nrows = 0; 63c4762a1bSJed Brown for (row=rstart; row<rend; row++) { 645f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRow(A,row,&ncols,&cols,&vals)); 65c4762a1bSJed Brown nnzA += ncols; 66c4762a1bSJed Brown norm = 0.0; 67c4762a1bSJed Brown for (j=0; j<ncols; j++) { 68c4762a1bSJed Brown val = PetscAbsScalar(vals[j]); 69c4762a1bSJed Brown if (norm < val) norm = norm; 70c4762a1bSJed Brown if (val > dtol) { 715f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(Asp,1,&row,1,&cols[j],&vals[j],INSERT_VALUES)); 72c4762a1bSJed Brown nnzAsp++; 73c4762a1bSJed Brown } 74c4762a1bSJed Brown } 75c4762a1bSJed Brown if (!norm) nrows++; 765f80ce2aSJacob Faibussowitsch CHKERRQ(MatRestoreRow(A,row,&ncols,&cols,&vals)); 77c4762a1bSJed Brown } 785f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(Asp,MAT_FINAL_ASSEMBLY)); 795f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(Asp,MAT_FINAL_ASSEMBLY)); 80c4762a1bSJed Brown 81c4762a1bSJed Brown percent=(PetscReal)nnzA*100/(m*n); 825f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF," [%d] Matrix A local size %d,%d; nnzA %d, %g percent; No. of zero rows: %d\n",rank,m,n,nnzA,percent,nrows)); 83c4762a1bSJed Brown percent=(PetscReal)nnzAsp*100/(m*n); 845f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF," [%d] Matrix Asp nnzAsp %d, %g percent\n",rank,nnzAsp,percent)); 85c4762a1bSJed Brown 86c4762a1bSJed Brown /* investigate matcoloring for Asp */ 87c4762a1bSJed Brown PetscBool Asp_coloring = PETSC_FALSE; 885f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-Asp_color",&Asp_coloring)); 89c4762a1bSJed Brown if (Asp_coloring) { 90c4762a1bSJed Brown MatColoring mc; 91c4762a1bSJed Brown ISColoring iscoloring; 92c4762a1bSJed Brown MatFDColoring matfdcoloring; 935f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD," Create coloring of Asp...\n")); 945f80ce2aSJacob Faibussowitsch CHKERRQ(MatColoringCreate(Asp,&mc)); 955f80ce2aSJacob Faibussowitsch CHKERRQ(MatColoringSetType(mc,MATCOLORINGSL)); 965f80ce2aSJacob Faibussowitsch CHKERRQ(MatColoringSetFromOptions(mc)); 975f80ce2aSJacob Faibussowitsch CHKERRQ(MatColoringApply(mc,&iscoloring)); 985f80ce2aSJacob Faibussowitsch CHKERRQ(MatColoringDestroy(&mc)); 995f80ce2aSJacob Faibussowitsch CHKERRQ(MatFDColoringCreate(Asp,iscoloring,&matfdcoloring)); 1005f80ce2aSJacob Faibussowitsch CHKERRQ(MatFDColoringSetFromOptions(matfdcoloring)); 1015f80ce2aSJacob Faibussowitsch CHKERRQ(MatFDColoringSetUp(Asp,iscoloring,matfdcoloring)); 1025f80ce2aSJacob Faibussowitsch /*CHKERRQ(MatFDColoringView(matfdcoloring,PETSC_VIEWER_STDOUT_WORLD));*/ 1035f80ce2aSJacob Faibussowitsch CHKERRQ(ISColoringDestroy(&iscoloring)); 1045f80ce2aSJacob Faibussowitsch CHKERRQ(MatFDColoringDestroy(&matfdcoloring)); 105c4762a1bSJed Brown } 106c4762a1bSJed Brown 107c4762a1bSJed Brown /* Write Asp in binary for study - see ~petsc/src/mat/tests/ex124.c */ 108c4762a1bSJed Brown PetscBool Asp_write = PETSC_FALSE; 1095f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-Asp_write",&Asp_write)); 110c4762a1bSJed Brown if (Asp_write) { 111c4762a1bSJed Brown PetscViewer viewer; 1125f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Write Asp into file Asp.dat ...\n")); 1135f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"Asp.dat",FILE_MODE_WRITE,&viewer)); 1145f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(Asp,viewer)); 1155f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&viewer)); 116c4762a1bSJed Brown } 117c4762a1bSJed Brown 1185f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 1195f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Asp)); 120*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 121*b122ec5aSJacob Faibussowitsch return 0; 122c4762a1bSJed Brown } 123