Actual source code: ex16.c
2: static char help[] = "Tests MatGetArray().\n\n";
4: #include petscmat.h
8: int main(int argc,char **args)
9: {
10: Mat A;
11: PetscInt i,j,m = 3,n = 2,rstart,rend;
12: PetscErrorCode ierr;
13: PetscScalar v,*array;
15: PetscInitialize(&argc,&args,(char *)0,help);
17: /*
18: Create a parallel dense matrix shared by all processors
19: */
20: MatCreateMPIDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,m,n,PETSC_NULL,&A);
21:
23: /*
24: Set values into the matrix
25: */
26: for (i=0; i<m; i++) {
27: for (j=0; j<n; j++) {
28: v = 9.0/(i+j+1); MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES);
29: }
30: }
31: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
32: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
34: /*
35: Print the matrix to the screen
36: */
37: MatView(A,PETSC_VIEWER_STDOUT_WORLD);
40: /*
41: Print the local portion of the matrix to the screen
42: */
43: MatGetArray(A,&array);
44: MatGetOwnershipRange(A,&rstart,&rend);
45: for (i=rstart; i<rend; i++) {
46: for (j=0; j<n; j++) {
47: PetscSynchronizedPrintf(PETSC_COMM_WORLD,"%6.4e ",PetscRealPart(array[j*(rend-rstart)+i-rstart]));
48: }
49: PetscSynchronizedPrintf(PETSC_COMM_WORLD,"\n");
50: }
51: PetscSynchronizedFlush(PETSC_COMM_WORLD);
52: MatRestoreArray(A,&array);
54: /*
55: Free the space used by the matrix
56: */
57: MatDestroy(A);
58: PetscFinalize();
59: return 0;
60: }