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