xref: /petsc/src/mat/tests/ex61.c (revision 40badf4fbc550ac1f60bd080eaff6de6d55b946d)
1 
2 static char help[] = "Tests MatSeq(B)AIJSetColumnIndices().\n\n";
3 
4 #include <petscmat.h>
5 
6 /*
7       Generate the following matrix:
8 
9          1 0 3
10          1 2 3
11          0 0 3
12 */
13 int main(int argc,char **args)
14 {
15   Mat            A;
16   PetscScalar    v;
17   PetscErrorCode ierr;
18   PetscInt       i,j,rowlens[] = {2,3,1},cols[] = {0,2,0,1,2,2};
19   PetscBool      flg;
20 
21   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
22   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-baij",&flg));
23   if (flg) {
24     CHKERRQ(MatCreateSeqBAIJ(PETSC_COMM_WORLD,1,3,3,0,rowlens,&A));
25     CHKERRQ(MatSeqBAIJSetColumnIndices(A,cols));
26   } else {
27     CHKERRQ(MatCreateSeqAIJ(PETSC_COMM_WORLD,3,3,0,rowlens,&A));
28     CHKERRQ(MatSeqAIJSetColumnIndices(A,cols));
29   }
30 
31   i    = 0; j = 0; v = 1.0;
32   CHKERRQ(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES));
33   i    = 0; j = 2; v = 3.0;
34   CHKERRQ(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES));
35 
36   i    = 1; j = 0; v = 1.0;
37   CHKERRQ(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES));
38   i    = 1; j = 1; v = 2.0;
39   CHKERRQ(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES));
40   i    = 1; j = 2; v = 3.0;
41   CHKERRQ(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES));
42 
43   i    = 2; j = 2; v = 3.0;
44   CHKERRQ(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES));
45 
46   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
47   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
48   CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
49 
50   CHKERRQ(MatDestroy(&A));
51   ierr = PetscFinalize();
52   return ierr;
53 }
54 
55 /*TEST
56 
57    test:
58 
59    test:
60       suffix: 2
61       args: -baij
62 
63 TEST*/
64