xref: /petsc/src/mat/tests/ex137.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests MatCreateMPISBAIJWithArrays().\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6*d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
7*d71ae5a4SJacob Faibussowitsch {
8c4762a1bSJed Brown   Mat         A;
9c4762a1bSJed Brown   PetscInt    m = 4, bs = 1, ii[5], jj[7];
10c4762a1bSJed Brown   PetscMPIInt size, rank;
11c4762a1bSJed Brown   PetscScalar aa[7];
12c4762a1bSJed Brown 
13327415f7SBarry Smith   PetscFunctionBeginUser;
149566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
159566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
169566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
17be096a46SBarry Smith   PetscCheck(size == 2, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for two processes");
18c4762a1bSJed Brown 
19dd400576SPatrick Sanan   if (rank == 0) {
209371c9d4SSatish Balay     ii[0] = 0;
219371c9d4SSatish Balay     ii[1] = 2;
229371c9d4SSatish Balay     ii[2] = 5;
239371c9d4SSatish Balay     ii[3] = 7;
249371c9d4SSatish Balay     ii[4] = 7;
259371c9d4SSatish Balay     jj[0] = 0;
269371c9d4SSatish Balay     jj[1] = 1;
279371c9d4SSatish Balay     jj[2] = 1;
289371c9d4SSatish Balay     jj[3] = 2;
299371c9d4SSatish Balay     jj[4] = 6;
309371c9d4SSatish Balay     jj[5] = 3;
319371c9d4SSatish Balay     jj[6] = 7;
329371c9d4SSatish Balay     aa[0] = 0;
339371c9d4SSatish Balay     aa[1] = 1;
349371c9d4SSatish Balay     aa[2] = 2;
359371c9d4SSatish Balay     aa[3] = 3;
369371c9d4SSatish Balay     aa[4] = 4;
379371c9d4SSatish Balay     aa[5] = 5;
389371c9d4SSatish Balay     aa[6] = 6;
39c4762a1bSJed Brown     /*  0 1
40c4762a1bSJed Brown           1  2       6
41c4762a1bSJed Brown              3          7 */
42c4762a1bSJed Brown   } else {
439371c9d4SSatish Balay     ii[0] = 0;
449371c9d4SSatish Balay     ii[1] = 2;
459371c9d4SSatish Balay     ii[2] = 4;
469371c9d4SSatish Balay     ii[3] = 6;
479371c9d4SSatish Balay     ii[4] = 7;
489371c9d4SSatish Balay     jj[0] = 4;
499371c9d4SSatish Balay     jj[1] = 5;
509371c9d4SSatish Balay     jj[2] = 5;
519371c9d4SSatish Balay     jj[3] = 7;
529371c9d4SSatish Balay     jj[4] = 6;
539371c9d4SSatish Balay     jj[5] = 7;
549371c9d4SSatish Balay     jj[6] = 7;
559371c9d4SSatish Balay     aa[0] = 8;
569371c9d4SSatish Balay     aa[1] = 9;
579371c9d4SSatish Balay     aa[2] = 10;
589371c9d4SSatish Balay     aa[3] = 11;
599371c9d4SSatish Balay     aa[4] = 12;
609371c9d4SSatish Balay     aa[5] = 13;
619371c9d4SSatish Balay     aa[6] = 14;
62c4762a1bSJed Brown     /*    4  5
63c4762a1bSJed Brown              5   7
64c4762a1bSJed Brown                6 7
65c4762a1bSJed Brown                  7 */
66c4762a1bSJed Brown   }
679566063dSJacob Faibussowitsch   PetscCall(MatCreateMPISBAIJWithArrays(PETSC_COMM_WORLD, bs, m, m, PETSC_DECIDE, PETSC_DECIDE, ii, jj, aa, &A));
689566063dSJacob Faibussowitsch   PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD));
699566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
709566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
71b122ec5aSJacob Faibussowitsch   return 0;
72c4762a1bSJed Brown }
73c4762a1bSJed Brown 
74c4762a1bSJed Brown /*TEST
75c4762a1bSJed Brown 
76c4762a1bSJed Brown    test:
77c4762a1bSJed Brown       nsize: 2
78c4762a1bSJed Brown 
79c4762a1bSJed Brown TEST*/
80