xref: /petsc/src/mat/tests/ex176.c (revision fdf6c4e30aafdbc795e4f76379caa977fd5cdf5a) !
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   Mat            A;
23   PetscInt       i[3][3] = {{0, 3, 6},{0, 3, 7},{0, 3, 5}};
24   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}};
25   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}};
26   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}};
27   MPI_Comm       comm;
28   PetscMPIInt    rank,size;
29 
30   PetscCall(PetscInitialize(&argc,&argv,NULL,help));
31   comm = PETSC_COMM_WORLD;
32   PetscCallMPI(MPI_Comm_rank(comm,&rank));
33   PetscCallMPI(MPI_Comm_size(comm,&size));
34   PetscCheck(size == 3,comm,PETSC_ERR_WRONG_MPI_SIZE,"You have to use three MPI processes to run this example ");
35   PetscCall(MatCreateMPIAIJWithArrays(comm,2,2,PETSC_DETERMINE,PETSC_DETERMINE,i[rank],j[rank],a[rank],&A));
36   PetscCall(MatView(A,NULL));
37   PetscCall(MatUpdateMPIAIJWithArrays(A,2,2,PETSC_DETERMINE,PETSC_DETERMINE,i[rank],j[rank],anew[rank]));
38   PetscCall(MatView(A,NULL));
39   PetscCall(MatUpdateMPIAIJWithArrays(A,2,2,PETSC_DETERMINE,PETSC_DETERMINE,i[rank],j[rank],a[rank]));
40   PetscCall(MatView(A,NULL));
41   PetscCall(MatUpdateMPIAIJWithArray(A,anew[rank]));
42   PetscCall(MatView(A,NULL));
43   PetscCall(MatDestroy(&A));
44   PetscCall(PetscFinalize());
45   return 0;
46 }
47 
48 /*TEST
49    test:
50      nsize: 3
51 
52 TEST*/
53