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