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