xref: /petsc/src/mat/tests/ex137.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests MatCreateMPISBAIJWithArrays().\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown int main(int argc,char **args)
7c4762a1bSJed Brown {
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 
13*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
145f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
155f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
162c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size != 2,PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for two processes");
17c4762a1bSJed Brown 
18dd400576SPatrick Sanan   if (rank == 0) {
19c4762a1bSJed Brown     ii[0] = 0; ii[1] = 2; ii[2] = 5; ii[3] = 7; ii[4] = 7;
20c4762a1bSJed Brown     jj[0] = 0; jj[1] = 1; jj[2] = 1; jj[3] = 2; jj[4] = 6; jj[5] = 3; jj[6] = 7;
21c4762a1bSJed Brown     aa[0] = 0; aa[1] = 1; aa[2] = 2; aa[3] = 3; aa[4] = 4; aa[5] = 5; aa[6] = 6;
22c4762a1bSJed Brown     /*  0 1
23c4762a1bSJed Brown           1  2       6
24c4762a1bSJed Brown              3          7 */
25c4762a1bSJed Brown   } else {
26c4762a1bSJed Brown     ii[0] = 0; ii[1] = 2; ii[2] = 4; ii[3] = 6; ii[4] = 7;
27c4762a1bSJed Brown     jj[0] = 4; jj[1] = 5; jj[2] = 5; jj[3] = 7; jj[4] = 6; jj[5] = 7; jj[6] = 7;
28c4762a1bSJed Brown     aa[0] = 8; aa[1] = 9; aa[2] = 10; aa[3] = 11; aa[4] = 12; aa[5] = 13; aa[6] = 14;
29c4762a1bSJed Brown     /*    4  5
30c4762a1bSJed Brown              5   7
31c4762a1bSJed Brown                6 7
32c4762a1bSJed Brown                  7 */
33c4762a1bSJed Brown   }
345f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateMPISBAIJWithArrays(PETSC_COMM_WORLD,bs,m,m,PETSC_DECIDE,PETSC_DECIDE,ii,jj,aa,&A));
355f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
365f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
37*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
38*b122ec5aSJacob Faibussowitsch   return 0;
39c4762a1bSJed Brown }
40c4762a1bSJed Brown 
41c4762a1bSJed Brown /*TEST
42c4762a1bSJed Brown 
43c4762a1bSJed Brown    test:
44c4762a1bSJed Brown       nsize: 2
45c4762a1bSJed Brown 
46c4762a1bSJed Brown TEST*/
47