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