1 static char help[] = "Passes a sparse matrix to MATLAB.\n\n"; 2 3 #include <petscmat.h> 4 5 int main(int argc, char **args) 6 { 7 PetscInt m = 4, n = 5, i, j, II, J; 8 PetscScalar one = 1.0, v; 9 Vec x; 10 Mat A; 11 12 PetscFunctionBeginUser; 13 PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 14 PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL)); 15 PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL)); 16 17 PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 18 PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, m * n, m * n)); 19 PetscCall(MatSetFromOptions(A)); 20 PetscCall(MatSetUp(A)); 21 22 for (i = 0; i < m; i++) { 23 for (j = 0; j < n; j++) { 24 v = -1.0; 25 II = j + n * i; 26 if (i > 0) { 27 J = II - n; 28 PetscCall(MatSetValues(A, 1, &II, 1, &J, &v, INSERT_VALUES)); 29 } 30 if (i < m - 1) { 31 J = II + n; 32 PetscCall(MatSetValues(A, 1, &II, 1, &J, &v, INSERT_VALUES)); 33 } 34 if (j > 0) { 35 J = II - 1; 36 PetscCall(MatSetValues(A, 1, &II, 1, &J, &v, INSERT_VALUES)); 37 } 38 if (j < n - 1) { 39 J = II + 1; 40 PetscCall(MatSetValues(A, 1, &II, 1, &J, &v, INSERT_VALUES)); 41 } 42 v = 4.0; 43 PetscCall(MatSetValues(A, 1, &II, 1, &II, &v, INSERT_VALUES)); 44 } 45 } 46 PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 47 PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 48 #if defined(PETSC_USE_SOCKET_VIEWER) 49 PetscCall(MatView(A, PETSC_VIEWER_SOCKET_WORLD)); 50 #endif 51 PetscCall(VecCreateSeq(PETSC_COMM_SELF, m, &x)); 52 PetscCall(VecSet(x, one)); 53 #if defined(PETSC_USE_SOCKET_VIEWER) 54 PetscCall(VecView(x, PETSC_VIEWER_SOCKET_WORLD)); 55 #endif 56 57 PetscCall(PetscSleep(30)); 58 59 PetscCall(VecDestroy(&x)); 60 PetscCall(MatDestroy(&A)); 61 PetscCall(PetscFinalize()); 62 return 0; 63 } 64 65 /*TEST 66 67 test: 68 TODO: cannot test with socket viewer 69 70 TEST*/ 71