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