Actual source code: ex16.c
petsc-3.6.4 2016-04-12
2: static char help[] = "Tests MatDenseGetArray() and MatView_SeqDense_Binary(), MatView_MPIDense_Binary().\n\n";
4: #include <petscmat.h>
5: #include <petscviewer.h>
9: int main(int argc,char **args)
10: {
11: Mat A;
12: PetscInt i,j,m = 3,n = 2,rstart,rend;
14: PetscScalar v,*array;
15: PetscViewer view;
17: PetscInitialize(&argc,&args,(char*)0,help);
19: /*
20: Create a parallel dense matrix shared by all processors
21: */
22: MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,m,n,NULL,&A);
24: /*
25: Set values into the matrix
26: */
27: for (i=0; i<m; i++) {
28: for (j=0; j<n; j++) {
29: v = 9.0/(i+j+1); MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES);
30: }
31: }
32: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
33: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
35: /*
36: Print the matrix to the screen
37: */
38: MatView(A,PETSC_VIEWER_STDOUT_WORLD);
41: /*
42: Print the local portion of the matrix to the screen
43: */
44: MatDenseGetArray(A,&array);
45: MatGetOwnershipRange(A,&rstart,&rend);
46: for (i=rstart; i<rend; i++) {
47: for (j=0; j<n; j++) {
48: PetscSynchronizedPrintf(PETSC_COMM_WORLD,"%6.4e ",(double)PetscRealPart(array[j*(rend-rstart)+i-rstart]));
49: }
50: PetscSynchronizedPrintf(PETSC_COMM_WORLD,"\n");
51: }
52: PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);
53: MatDenseRestoreArray(A,&array);
55: /*
56: Store the binary matrix to a file
57: */
58: PetscViewerBinaryOpen(PETSC_COMM_WORLD, "matrix.dat", FILE_MODE_WRITE, &view);
59: PetscViewerSetFormat(view,PETSC_VIEWER_NATIVE);
60: MatView(A,view);
61: PetscViewerDestroy(&view);
62: MatDestroy(&A);
64: /*
65: Now reload the matrix and view it
66: */
67: PetscViewerBinaryOpen(PETSC_COMM_WORLD,"matrix.dat",FILE_MODE_READ,&view);
68: MatCreate(PETSC_COMM_WORLD,&A);
69: MatSetType(A,MATMPIDENSE);
70: MatLoad(A,view);
71: PetscViewerDestroy(&view);
72: MatView(A,PETSC_VIEWER_STDOUT_WORLD);
73: MatDestroy(&A);
75: PetscFinalize();
76: return 0;
77: }