Actual source code: ex16.c
petsc-3.10.5 2019-03-28
2: static char help[] = "Tests MatDenseGetArray() and MatView_SeqDense_Binary(), MatView_MPIDense_Binary().\n\n";
4: #include <petscmat.h>
5: #include <petscviewer.h>
7: int main(int argc,char **args)
8: {
9: Mat A;
10: PetscInt i,j,m = 3,n = 2,rstart,rend;
12: PetscScalar v,*array;
13: PetscViewer view;
15: PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
16: /*
17: Create a parallel dense matrix shared by all processors
18: */
19: MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,m,n,NULL,&A);
21: /*
22: Set values into the matrix
23: */
24: for (i=0; i<m; i++) {
25: for (j=0; j<n; j++) {
26: v = 9.0/(i+j+1); MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES);
27: }
28: }
29: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
30: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
32: /*
33: Print the matrix to the screen
34: */
35: MatView(A,PETSC_VIEWER_STDOUT_WORLD);
38: /*
39: Print the local portion of the matrix to the screen
40: */
41: MatDenseGetArray(A,&array);
42: MatGetOwnershipRange(A,&rstart,&rend);
43: for (i=rstart; i<rend; i++) {
44: for (j=0; j<n; j++) {
45: PetscSynchronizedPrintf(PETSC_COMM_WORLD,"%6.4e ",(double)PetscRealPart(array[j*(rend-rstart)+i-rstart]));
46: }
47: PetscSynchronizedPrintf(PETSC_COMM_WORLD,"\n");
48: }
49: PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);
50: MatDenseRestoreArray(A,&array);
52: /*
53: Store the binary matrix to a file
54: */
55: PetscViewerBinaryOpen(PETSC_COMM_WORLD, "matrix.dat", FILE_MODE_WRITE, &view);
56: PetscViewerPushFormat(view,PETSC_VIEWER_NATIVE);
57: MatView(A,view);
58: PetscViewerPopFormat(view);
59: PetscViewerDestroy(&view);
60: MatDestroy(&A);
62: /*
63: Now reload the matrix and view it
64: */
65: PetscViewerBinaryOpen(PETSC_COMM_WORLD,"matrix.dat",FILE_MODE_READ,&view);
66: MatCreate(PETSC_COMM_WORLD,&A);
67: MatSetType(A,MATMPIDENSE);
68: MatLoad(A,view);
69: PetscViewerDestroy(&view);
70: MatView(A,PETSC_VIEWER_STDOUT_WORLD);
73: PetscMalloc1((rend-rstart)*n,&array);
74: for (i=0; i<(rend-rstart)*n; i++) array[i] = 1.;
75: MatDensePlaceArray(A,array);
76: MatView(A,PETSC_VIEWER_STDOUT_WORLD);
77: MatDenseResetArray(A);
78: PetscFree(array);
79: MatView(A,PETSC_VIEWER_STDOUT_WORLD);
80: MatDestroy(&A);
82: PetscFinalize();
83: return ierr;
84: }
87: /*TEST
89: test:
90: nsize: 2
92: TEST*/