Actual source code: ex16.c
petsc-3.8.4 2018-03-24
2: static char help[] = "Reads a matrix from PETSc binary file. Use for view or investigating matrix data structure. \n\n";
3: /*
4: Example:
5: ./ex16 -f <matrix file> -a_mat_view draw -draw_pause -1
6: ./ex16 -f <matrix file> -a_mat_view ascii::ascii_info
7: */
9: #include <petscmat.h>
10: int main(int argc,char **args)
11: {
12: Mat A,Asp;
13: PetscViewer fd; /* viewer */
14: char file[PETSC_MAX_PATH_LEN]; /* input file name */
15: PetscErrorCode ierr;
16: PetscInt m,n,rstart,rend;
17: PetscBool flg;
18: PetscInt row,ncols,j,nrows,nnzA=0,nnzAsp=0;
19: const PetscInt *cols;
20: const PetscScalar *vals;
21: PetscReal norm,percent,val,dtol=1.e-16;
22: PetscMPIInt rank;
23: MatInfo matinfo;
24: PetscInt Dnnz,Onnz;
27: PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
28: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
30: /* Determine files from which we read the linear systems. */
31: PetscOptionsGetString(NULL,NULL,"-f",file,PETSC_MAX_PATH_LEN,&flg);
32: if (!flg) SETERRQ(PETSC_COMM_WORLD,1,"Must indicate binary file with the -f option");
34: /* Open binary file. Note that we use FILE_MODE_READ to indicate
35: reading from this file. */
36: PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd);
38: /* Load the matrix; then destroy the viewer. */
39: MatCreate(PETSC_COMM_WORLD,&A);
40: MatSetOptionsPrefix(A,"a_");
41: MatSetFromOptions(A);
42: MatLoad(A,fd);
43: PetscViewerDestroy(&fd);
44: MatGetSize(A,&m,&n);
45: MatGetInfo(A,MAT_LOCAL,&matinfo);
46: /*printf("matinfo.nz_used %g\n",matinfo.nz_used);*/
48: /* Get a sparse matrix Asp by dumping zero entries of A */
49: MatCreate(PETSC_COMM_WORLD,&Asp);
50: MatSetSizes(Asp,m,n,PETSC_DECIDE,PETSC_DECIDE);
51: MatSetOptionsPrefix(Asp,"asp_");
52: MatSetFromOptions(Asp);
53: Dnnz = (PetscInt)matinfo.nz_used/m + 1;
54: Onnz = Dnnz/2;
55: printf("Dnnz %d %d\n",Dnnz,Onnz);
56: MatSeqAIJSetPreallocation(Asp,Dnnz,NULL);
57: MatMPIAIJSetPreallocation(Asp,Dnnz,NULL,Onnz,NULL);
58: /* The allocation above is approximate so we must set this option to be permissive.
59: * Real code should preallocate exactly. */
60: MatSetOption(Asp,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);
62: /* Check zero rows */
63: MatGetOwnershipRange(A,&rstart,&rend);
64: nrows = 0;
65: for (row=rstart; row<rend; row++) {
66: MatGetRow(A,row,&ncols,&cols,&vals);
67: nnzA += ncols;
68: norm = 0.0;
69: for (j=0; j<ncols; j++) {
70: val = PetscAbsScalar(vals[j]);
71: if (norm < val) norm = norm;
72: if (val > dtol) {
73: MatSetValues(Asp,1,&row,1,&cols[j],&vals[j],INSERT_VALUES);
74: nnzAsp++;
75: }
76: }
77: if (!norm) nrows++;
78: MatRestoreRow(A,row,&ncols,&cols,&vals);
79: }
80: MatAssemblyBegin(Asp,MAT_FINAL_ASSEMBLY);
81: MatAssemblyEnd(Asp,MAT_FINAL_ASSEMBLY);
83: percent=(PetscReal)nnzA*100/(m*n);
84: 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);
85: percent=(PetscReal)nnzAsp*100/(m*n);
86: PetscPrintf(PETSC_COMM_SELF," [%d] Matrix Asp nnzAsp %d, %g percent\n",rank,nnzAsp,percent);
88: /* investigate matcoloring for Asp */
89: PetscBool Asp_coloring = PETSC_FALSE;
90: PetscOptionsHasName(NULL,NULL,"-Asp_color",&Asp_coloring);
91: if (Asp_coloring) {
92: MatColoring mc;
93: ISColoring iscoloring;
94: MatFDColoring matfdcoloring;
95: PetscPrintf(PETSC_COMM_WORLD," Create coloring of Asp...\n");
96: MatColoringCreate(Asp,&mc);
97: MatColoringSetType(mc,MATCOLORINGSL);
98: MatColoringSetFromOptions(mc);
99: MatColoringApply(mc,&iscoloring);
100: MatColoringDestroy(&mc);
101: MatFDColoringCreate(Asp,iscoloring,&matfdcoloring);
102: MatFDColoringSetFromOptions(matfdcoloring);
103: MatFDColoringSetUp(Asp,iscoloring,matfdcoloring);
104: /*MatFDColoringView(matfdcoloring,PETSC_VIEWER_STDOUT_WORLD);*/
105: ISColoringDestroy(&iscoloring);
106: MatFDColoringDestroy(&matfdcoloring);
107: }
109: /* Write Asp in binary for study - see ~petsc/src/mat/examples/tests/ex124.c */
110: PetscBool Asp_write = PETSC_FALSE;
111: PetscOptionsHasName(NULL,NULL,"-Asp_write",&Asp_write);
112: if (Asp_write) {
113: PetscViewer viewer;
114: PetscPrintf(PETSC_COMM_SELF,"Write Asp into file Asp.dat ...\n");
115: PetscViewerBinaryOpen(PETSC_COMM_WORLD,"Asp.dat",FILE_MODE_WRITE,&viewer);
116: MatView(Asp,viewer);
117: PetscViewerDestroy(&viewer);
118: }
120: MatDestroy(&A);
121: MatDestroy(&Asp);
122: PetscFinalize();
123: return ierr;
124: }