xref: /petsc/src/mat/tests/ex176.c (revision a0bcfa1fc66e6ceb348a18954676b654f6a1dad4)
1c4762a1bSJed Brown 
2*a0bcfa1fSJunchao Zhang static char help[] = "Tests MatCreateMPIAIJWithArrays() abd MatUpdateMPIAIJWithArray()\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown /*
7c4762a1bSJed Brown  * This is an extremely simple example to test MatUpdateMPIAIJWithArrays()
8c4762a1bSJed Brown  *
9c4762a1bSJed Brown  * A =
10c4762a1bSJed Brown 
11c4762a1bSJed Brown    1    2   0   3  0  0
12c4762a1bSJed Brown    0    4   5   0  0  6
13c4762a1bSJed Brown    7    0   8   0  9  0
14c4762a1bSJed Brown    0   10  11  12  0  13
15c4762a1bSJed Brown    0   14  15   0  0  16
16c4762a1bSJed Brown   17    0   0   0  0  18
17c4762a1bSJed Brown  *
18c4762a1bSJed Brown  * */
19c4762a1bSJed Brown 
20d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
21d71ae5a4SJacob Faibussowitsch {
22*a0bcfa1fSJunchao Zhang   Mat      A, B;
239371c9d4SSatish Balay   PetscInt i[3][3] = {
249371c9d4SSatish Balay     {0, 3, 6},
259371c9d4SSatish Balay     {0, 3, 7},
269371c9d4SSatish Balay     {0, 3, 5}
279371c9d4SSatish Balay   };
289371c9d4SSatish Balay   PetscInt j[3][7] = {
299371c9d4SSatish Balay     {0, 1, 3, 1, 2, 5,  -1},
309371c9d4SSatish Balay     {0, 2, 4, 1, 2, 3,  5 },
319371c9d4SSatish Balay     {1, 2, 5, 0, 5, -1, -1}
329371c9d4SSatish Balay   };
339371c9d4SSatish Balay   PetscScalar a[3][7] = {
349371c9d4SSatish Balay     {1,  2,  3,  4,  5,  6,  -1},
359371c9d4SSatish Balay     {7,  8,  9,  10, 11, 12, 13},
369371c9d4SSatish Balay     {14, 15, 16, 17, 18, -1, -1}
379371c9d4SSatish Balay   };
389371c9d4SSatish Balay   PetscScalar anew[3][7] = {
399371c9d4SSatish Balay     {10,  20,  30,  40,  50,  60,  -1 },
409371c9d4SSatish Balay     {70,  80,  90,  100, 110, 120, 130},
419371c9d4SSatish Balay     {140, 150, 160, 170, 180, -1,  -1 }
429371c9d4SSatish Balay   };
43c4762a1bSJed Brown   MPI_Comm    comm;
44*a0bcfa1fSJunchao Zhang   PetscMPIInt rank;
45*a0bcfa1fSJunchao Zhang   PetscBool   equal = PETSC_FALSE;
46c4762a1bSJed Brown 
47327415f7SBarry Smith   PetscFunctionBeginUser;
489566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
49c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
509566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
51*a0bcfa1fSJunchao Zhang   PetscCall(MatCreateMPIAIJWithArrays(comm, 2, PETSC_DETERMINE, PETSC_DETERMINE, 6, i[rank], j[rank], a[rank], &B));
52*a0bcfa1fSJunchao Zhang 
53*a0bcfa1fSJunchao Zhang   PetscCall(MatCreateMPIAIJWithArrays(comm, 2, PETSC_DETERMINE, PETSC_DETERMINE, 6, i[rank], j[rank], a[rank], &A));
54*a0bcfa1fSJunchao Zhang   PetscCall(MatSetFromOptions(A)); /* might change A's type */
55*a0bcfa1fSJunchao Zhang 
56*a0bcfa1fSJunchao Zhang   PetscCall(MatEqual(A, B, &equal));
57*a0bcfa1fSJunchao Zhang   PetscCheck(equal, comm, PETSC_ERR_PLIB, "wrong results");
58*a0bcfa1fSJunchao Zhang 
596a3d2595SBarry Smith   PetscCall(MatUpdateMPIAIJWithArray(A, anew[rank]));
60*a0bcfa1fSJunchao Zhang   PetscCall(MatUpdateMPIAIJWithArray(B, anew[rank]));
61*a0bcfa1fSJunchao Zhang   PetscCall(MatEqual(A, B, &equal));
62*a0bcfa1fSJunchao Zhang   PetscCheck(equal, comm, PETSC_ERR_PLIB, "wrong results");
63*a0bcfa1fSJunchao Zhang 
64*a0bcfa1fSJunchao Zhang   PetscCall(MatUpdateMPIAIJWithArray(A, a[rank]));
65*a0bcfa1fSJunchao Zhang   PetscCall(MatUpdateMPIAIJWithArray(B, a[rank]));
66*a0bcfa1fSJunchao Zhang   PetscCall(MatEqual(A, B, &equal));
67*a0bcfa1fSJunchao Zhang   PetscCheck(equal, comm, PETSC_ERR_PLIB, "wrong results");
68*a0bcfa1fSJunchao Zhang 
699566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
70*a0bcfa1fSJunchao Zhang   PetscCall(MatDestroy(&B));
719566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
72b122ec5aSJacob Faibussowitsch   return 0;
73c4762a1bSJed Brown }
74c4762a1bSJed Brown 
75c4762a1bSJed Brown /*TEST
76*a0bcfa1fSJunchao Zhang    testset:
77*a0bcfa1fSJunchao Zhang      nsize: {{1 3}}
78*a0bcfa1fSJunchao Zhang      output_file: output/empty.out
79c4762a1bSJed Brown 
80*a0bcfa1fSJunchao Zhang      test:
81*a0bcfa1fSJunchao Zhang        suffix: aij
82*a0bcfa1fSJunchao Zhang 
83*a0bcfa1fSJunchao Zhang      test:
84*a0bcfa1fSJunchao Zhang        requires: cuda
85*a0bcfa1fSJunchao Zhang        suffix: cuda
86*a0bcfa1fSJunchao Zhang        # since the matrices are created with MatCreateMPIxxx(), users are allowed to pass 'mpiaijcusparse' even with one rank
87*a0bcfa1fSJunchao Zhang        args: -mat_type {{aijcusparse mpiaijcusparse}}
88*a0bcfa1fSJunchao Zhang 
89*a0bcfa1fSJunchao Zhang      test:
90*a0bcfa1fSJunchao Zhang        requires: kokkos_kernels
91*a0bcfa1fSJunchao Zhang        suffix: kok
92*a0bcfa1fSJunchao Zhang        args: -mat_type aijkokkos
93c4762a1bSJed Brown TEST*/
94