xref: /petsc/src/mat/tests/ex176.c (revision 030f984af8d8bb4c203755d35bded3c05b3d83ce)
1 
2 static char help[] ="Tests MatCreateMPIAIJWithArrays() abd MatUpdateMPIAIJWithArrays()\n";
3 
4 #include <petscmat.h>
5 
6 /*
7  * This is an extremely simple example to test MatUpdateMPIAIJWithArrays()
8  *
9  * A =
10 
11    1    2   0   3  0  0
12    0    4   5   0  0  6
13    7    0   8   0  9  0
14    0   10  11  12  0  13
15    0   14  15   0  0  16
16   17    0   0   0  0  18
17  *
18  * */
19 
20 int main(int argc,char **argv)
21 {
22   PetscErrorCode ierr;
23   Mat            A;
24   PetscInt       i[3][3] = {{0, 3, 6},{0, 3, 7},{0, 3, 5}};
25   PetscInt       j[3][7] = {{0, 1, 3, 1, 2, 5, -1},{0, 2, 4, 1, 2, 3, 5},{1, 2, 5, 0, 5, -1, -1}};
26   PetscScalar    a[3][7] = {{1, 2, 3, 4, 5, 6, -1}, {7, 8, 9, 10, 11, 12, 13},{14, 15, 16, 17, 18, -1, -1}};
27   PetscScalar    anew[3][7] = {{10, 20, 30, 40, 50, 60, -1}, {70, 80, 90, 100, 110, 120, 130},{140, 150, 160, 170, 180, -1, -1}};
28   MPI_Comm       comm;
29   PetscMPIInt    rank,size;
30 
31   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
32   comm = PETSC_COMM_WORLD;
33   ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr);
34   ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr);
35   if (size != 3) SETERRQ(comm,PETSC_ERR_ARG_INCOMP,"You have to use three MPI processes to run this example \n");
36   ierr = MatCreateMPIAIJWithArrays(comm,2,2,PETSC_DETERMINE,PETSC_DETERMINE,i[rank],j[rank],a[rank],&A);CHKERRQ(ierr);
37   ierr = MatView(A,NULL);CHKERRQ(ierr);
38   ierr = MatUpdateMPIAIJWithArrays(A,2,2,PETSC_DETERMINE,PETSC_DETERMINE,i[rank],j[rank],anew[rank]);CHKERRQ(ierr);
39   ierr = MatView(A,NULL);CHKERRQ(ierr);
40   ierr = MatUpdateMPIAIJWithArrays(A,2,2,PETSC_DETERMINE,PETSC_DETERMINE,i[rank],j[rank],a[rank]);CHKERRQ(ierr);
41   ierr = MatView(A,NULL);CHKERRQ(ierr);
42   ierr = MatDestroy(&A);CHKERRQ(ierr);
43   ierr = PetscFinalize();
44   return ierr;
45 }
46 
47 /*TEST
48    test:
49      nsize: 3
50 
51 TEST*/
52