xref: /petsc/src/mat/tests/ex176.c (revision df4cd43f92eaa320656440c40edb1046daee8f75)
1 
2 static char help[] = "Tests MatCreateMPIAIJWithArrays() abd MatUpdateMPIAIJWithArray()\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 {
22   Mat      A, B;
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;
45   PetscBool   equal = PETSC_FALSE;
46 
47   PetscFunctionBeginUser;
48   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
49   comm = PETSC_COMM_WORLD;
50   PetscCallMPI(MPI_Comm_rank(comm, &rank));
51   PetscCall(MatCreateMPIAIJWithArrays(comm, 2, PETSC_DETERMINE, PETSC_DETERMINE, 6, i[rank], j[rank], a[rank], &B));
52 
53   PetscCall(MatCreateMPIAIJWithArrays(comm, 2, PETSC_DETERMINE, PETSC_DETERMINE, 6, i[rank], j[rank], a[rank], &A));
54   PetscCall(MatSetFromOptions(A)); /* might change A's type */
55 
56   PetscCall(MatEqual(A, B, &equal));
57   PetscCheck(equal, comm, PETSC_ERR_PLIB, "wrong results");
58 
59   PetscCall(MatUpdateMPIAIJWithArray(A, anew[rank]));
60   PetscCall(MatUpdateMPIAIJWithArray(B, anew[rank]));
61   PetscCall(MatEqual(A, B, &equal));
62   PetscCheck(equal, comm, PETSC_ERR_PLIB, "wrong results");
63 
64   PetscCall(MatUpdateMPIAIJWithArray(A, a[rank]));
65   PetscCall(MatUpdateMPIAIJWithArray(B, a[rank]));
66   PetscCall(MatEqual(A, B, &equal));
67   PetscCheck(equal, comm, PETSC_ERR_PLIB, "wrong results");
68 
69   PetscCall(MatDestroy(&A));
70   PetscCall(MatDestroy(&B));
71   PetscCall(PetscFinalize());
72   return 0;
73 }
74 
75 /*TEST
76    testset:
77      nsize: {{1 3}}
78      output_file: output/empty.out
79 
80      test:
81        suffix: aij
82 
83      test:
84        requires: cuda
85        suffix: cuda
86        # since the matrices are created with MatCreateMPIxxx(), users are allowed to pass 'mpiaijcusparse' even with one rank
87        args: -mat_type {{aijcusparse mpiaijcusparse}}
88 
89      test:
90        requires: kokkos_kernels
91        suffix: kok
92        args: -mat_type aijkokkos
93 TEST*/
94