Actual source code: ex71.c
petsc-3.6.1 2015-08-06
2: static char help[] = "Passes a sparse matrix to MATLAB.\n\n";
4: #include <petscmat.h>
8: int main(int argc,char **args)
9: {
11: PetscInt m = 4,n = 5,i,j,II,J;
12: PetscScalar one = 1.0,v;
13: Vec x;
14: Mat A;
16: PetscInitialize(&argc,&args,(char*)0,help);
17: PetscOptionsGetInt(NULL,"-m",&m,NULL);
18: PetscOptionsGetInt(NULL,"-n",&n,NULL);
20: MatCreate(PETSC_COMM_WORLD,&A);
21: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);
22: MatSetFromOptions(A);
24: for (i=0; i<m; i++) {
25: for (j=0; j<n; j++) {
26: v = -1.0; II = j + n*i;
27: if (i>0) {J = II - n; MatSetValues(A,1,&II,1,&J,&v,INSERT_VALUES);}
28: if (i<m-1) {J = II + n; MatSetValues(A,1,&II,1,&J,&v,INSERT_VALUES);}
29: if (j>0) {J = II - 1; MatSetValues(A,1,&II,1,&J,&v,INSERT_VALUES);}
30: if (j<n-1) {J = II + 1; MatSetValues(A,1,&II,1,&J,&v,INSERT_VALUES);}
31: v = 4.0; MatSetValues(A,1,&II,1,&II,&v,INSERT_VALUES);
32: }
33: }
34: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
35: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
36: #if defined(PETSC_USE_SOCKET_VIEWER)
37: MatView(A,PETSC_VIEWER_SOCKET_WORLD);
38: #endif
39: VecCreateSeq(PETSC_COMM_SELF,m,&x);
40: VecSet(x,one);
41: #if defined(PETSC_USE_SOCKET_VIEWER)
42: VecView(x,PETSC_VIEWER_SOCKET_WORLD);
43: #endif
45: PetscSleep(30);
47: VecDestroy(&x);
48: MatDestroy(&A);
49: PetscFinalize();
50: return 0;
51: }