xref: /petsc/src/mat/tests/ex137.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests MatCreateMPISBAIJWithArrays().\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6*9371c9d4SSatish Balay int main(int argc, char **args) {
7c4762a1bSJed Brown   Mat         A;
8c4762a1bSJed Brown   PetscInt    m = 4, bs = 1, ii[5], jj[7];
9c4762a1bSJed Brown   PetscMPIInt size, rank;
10c4762a1bSJed Brown   PetscScalar aa[7];
11c4762a1bSJed Brown 
12327415f7SBarry Smith   PetscFunctionBeginUser;
139566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
149566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
159566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
16be096a46SBarry Smith   PetscCheck(size == 2, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for two processes");
17c4762a1bSJed Brown 
18dd400576SPatrick Sanan   if (rank == 0) {
19*9371c9d4SSatish Balay     ii[0] = 0;
20*9371c9d4SSatish Balay     ii[1] = 2;
21*9371c9d4SSatish Balay     ii[2] = 5;
22*9371c9d4SSatish Balay     ii[3] = 7;
23*9371c9d4SSatish Balay     ii[4] = 7;
24*9371c9d4SSatish Balay     jj[0] = 0;
25*9371c9d4SSatish Balay     jj[1] = 1;
26*9371c9d4SSatish Balay     jj[2] = 1;
27*9371c9d4SSatish Balay     jj[3] = 2;
28*9371c9d4SSatish Balay     jj[4] = 6;
29*9371c9d4SSatish Balay     jj[5] = 3;
30*9371c9d4SSatish Balay     jj[6] = 7;
31*9371c9d4SSatish Balay     aa[0] = 0;
32*9371c9d4SSatish Balay     aa[1] = 1;
33*9371c9d4SSatish Balay     aa[2] = 2;
34*9371c9d4SSatish Balay     aa[3] = 3;
35*9371c9d4SSatish Balay     aa[4] = 4;
36*9371c9d4SSatish Balay     aa[5] = 5;
37*9371c9d4SSatish Balay     aa[6] = 6;
38c4762a1bSJed Brown     /*  0 1
39c4762a1bSJed Brown           1  2       6
40c4762a1bSJed Brown              3          7 */
41c4762a1bSJed Brown   } else {
42*9371c9d4SSatish Balay     ii[0] = 0;
43*9371c9d4SSatish Balay     ii[1] = 2;
44*9371c9d4SSatish Balay     ii[2] = 4;
45*9371c9d4SSatish Balay     ii[3] = 6;
46*9371c9d4SSatish Balay     ii[4] = 7;
47*9371c9d4SSatish Balay     jj[0] = 4;
48*9371c9d4SSatish Balay     jj[1] = 5;
49*9371c9d4SSatish Balay     jj[2] = 5;
50*9371c9d4SSatish Balay     jj[3] = 7;
51*9371c9d4SSatish Balay     jj[4] = 6;
52*9371c9d4SSatish Balay     jj[5] = 7;
53*9371c9d4SSatish Balay     jj[6] = 7;
54*9371c9d4SSatish Balay     aa[0] = 8;
55*9371c9d4SSatish Balay     aa[1] = 9;
56*9371c9d4SSatish Balay     aa[2] = 10;
57*9371c9d4SSatish Balay     aa[3] = 11;
58*9371c9d4SSatish Balay     aa[4] = 12;
59*9371c9d4SSatish Balay     aa[5] = 13;
60*9371c9d4SSatish Balay     aa[6] = 14;
61c4762a1bSJed Brown     /*    4  5
62c4762a1bSJed Brown              5   7
63c4762a1bSJed Brown                6 7
64c4762a1bSJed Brown                  7 */
65c4762a1bSJed Brown   }
669566063dSJacob Faibussowitsch   PetscCall(MatCreateMPISBAIJWithArrays(PETSC_COMM_WORLD, bs, m, m, PETSC_DECIDE, PETSC_DECIDE, ii, jj, aa, &A));
679566063dSJacob Faibussowitsch   PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD));
689566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
699566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
70b122ec5aSJacob Faibussowitsch   return 0;
71c4762a1bSJed Brown }
72c4762a1bSJed Brown 
73c4762a1bSJed Brown /*TEST
74c4762a1bSJed Brown 
75c4762a1bSJed Brown    test:
76c4762a1bSJed Brown       nsize: 2
77c4762a1bSJed Brown 
78c4762a1bSJed Brown TEST*/
79