xref: /petsc/src/mat/tests/ex137.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
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   PetscErrorCode ierr;
10c4762a1bSJed Brown   PetscInt       m = 4, bs = 1,ii[5],jj[7];
11c4762a1bSJed Brown   PetscMPIInt    size,rank;
12c4762a1bSJed Brown   PetscScalar    aa[7];
13c4762a1bSJed Brown 
14c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
15ffc4695bSBarry Smith   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
16ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr);
17*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size != 2,PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for two processes");
18c4762a1bSJed Brown 
19dd400576SPatrick Sanan   if (rank == 0) {
20c4762a1bSJed Brown     ii[0] = 0; ii[1] = 2; ii[2] = 5; ii[3] = 7; ii[4] = 7;
21c4762a1bSJed Brown     jj[0] = 0; jj[1] = 1; jj[2] = 1; jj[3] = 2; jj[4] = 6; jj[5] = 3; jj[6] = 7;
22c4762a1bSJed Brown     aa[0] = 0; aa[1] = 1; aa[2] = 2; aa[3] = 3; aa[4] = 4; aa[5] = 5; aa[6] = 6;
23c4762a1bSJed Brown     /*  0 1
24c4762a1bSJed Brown           1  2       6
25c4762a1bSJed Brown              3          7 */
26c4762a1bSJed Brown   } else {
27c4762a1bSJed Brown     ii[0] = 0; ii[1] = 2; ii[2] = 4; ii[3] = 6; ii[4] = 7;
28c4762a1bSJed Brown     jj[0] = 4; jj[1] = 5; jj[2] = 5; jj[3] = 7; jj[4] = 6; jj[5] = 7; jj[6] = 7;
29c4762a1bSJed Brown     aa[0] = 8; aa[1] = 9; aa[2] = 10; aa[3] = 11; aa[4] = 12; aa[5] = 13; aa[6] = 14;
30c4762a1bSJed Brown     /*    4  5
31c4762a1bSJed Brown              5   7
32c4762a1bSJed Brown                6 7
33c4762a1bSJed Brown                  7 */
34c4762a1bSJed Brown   }
35c4762a1bSJed Brown   ierr = MatCreateMPISBAIJWithArrays(PETSC_COMM_WORLD,bs,m,m,PETSC_DECIDE,PETSC_DECIDE,ii,jj,aa,&A);CHKERRQ(ierr);
36c4762a1bSJed Brown   ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
37c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
38c4762a1bSJed Brown   ierr = PetscFinalize();
39c4762a1bSJed Brown   return ierr;
40c4762a1bSJed Brown }
41c4762a1bSJed Brown 
42c4762a1bSJed Brown /*TEST
43c4762a1bSJed Brown 
44c4762a1bSJed Brown    test:
45c4762a1bSJed Brown       nsize: 2
46c4762a1bSJed Brown 
47c4762a1bSJed Brown TEST*/
48