Actual source code: ex176.c
2: static char help[] = "Tests MatCreateMPIAIJWithArrays() abd MatUpdateMPIAIJWithArrays()\n";
4: #include <petscmat.h>
6: /*
7: * This is an extremely simple example to test MatUpdateMPIAIJWithArrays()
8: *
9: * A =
11: 1 2 0 3 0 0
12: 0 4 5 0 0 6
13: 7 0 8 0 9 0
14: 0 10 11 12 0 13
15: 0 14 15 0 0 16
16: 17 0 0 0 0 18
17: *
18: * */
20: int main(int argc, char **argv)
21: {
22: Mat A;
23: PetscInt i[3][3] = {
24: {0, 3, 6},
25: {0, 3, 7},
26: {0, 3, 5}
27: };
28: PetscInt j[3][7] = {
29: {0, 1, 3, 1, 2, 5, -1},
30: {0, 2, 4, 1, 2, 3, 5 },
31: {1, 2, 5, 0, 5, -1, -1}
32: };
33: PetscScalar a[3][7] = {
34: {1, 2, 3, 4, 5, 6, -1},
35: {7, 8, 9, 10, 11, 12, 13},
36: {14, 15, 16, 17, 18, -1, -1}
37: };
38: PetscScalar anew[3][7] = {
39: {10, 20, 30, 40, 50, 60, -1 },
40: {70, 80, 90, 100, 110, 120, 130},
41: {140, 150, 160, 170, 180, -1, -1 }
42: };
43: MPI_Comm comm;
44: PetscMPIInt rank, size;
47: PetscInitialize(&argc, &argv, NULL, help);
48: comm = PETSC_COMM_WORLD;
49: MPI_Comm_rank(comm, &rank);
50: MPI_Comm_size(comm, &size);
52: MatCreateMPIAIJWithArrays(comm, 2, 2, PETSC_DETERMINE, PETSC_DETERMINE, i[rank], j[rank], a[rank], &A);
53: MatView(A, NULL);
54: MatUpdateMPIAIJWithArrays(A, 2, 2, PETSC_DETERMINE, PETSC_DETERMINE, i[rank], j[rank], anew[rank]);
55: MatView(A, NULL);
56: MatUpdateMPIAIJWithArrays(A, 2, 2, PETSC_DETERMINE, PETSC_DETERMINE, i[rank], j[rank], a[rank]);
57: MatView(A, NULL);
58: MatUpdateMPIAIJWithArray(A, anew[rank]);
59: MatView(A, NULL);
60: MatDestroy(&A);
61: PetscFinalize();
62: return 0;
63: }
65: /*TEST
66: test:
67: nsize: 3
69: TEST*/