xref: /petsc/src/mat/tests/ex61.c (revision ffa8c5705e8ab2cf85ee1d14dbe507a6e2eb5283)
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   PetscInt       i,j,rowlens[] = {2,3,1},cols[] = {0,2,0,1,2,2};
18   PetscBool      flg;
19 
20   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
21   PetscCall(PetscOptionsHasName(NULL,NULL,"-baij",&flg));
22   if (flg) {
23     PetscCall(MatCreateSeqBAIJ(PETSC_COMM_WORLD,1,3,3,0,rowlens,&A));
24     PetscCall(MatSeqBAIJSetColumnIndices(A,cols));
25   } else {
26     PetscCall(MatCreateSeqAIJ(PETSC_COMM_WORLD,3,3,0,rowlens,&A));
27     PetscCall(MatSeqAIJSetColumnIndices(A,cols));
28   }
29 
30   i    = 0; j = 0; v = 1.0;
31   PetscCall(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES));
32   i    = 0; j = 2; v = 3.0;
33   PetscCall(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES));
34 
35   i    = 1; j = 0; v = 1.0;
36   PetscCall(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES));
37   i    = 1; j = 1; v = 2.0;
38   PetscCall(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES));
39   i    = 1; j = 2; v = 3.0;
40   PetscCall(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES));
41 
42   i    = 2; j = 2; v = 3.0;
43   PetscCall(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES));
44 
45   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
46   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
47   PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
48 
49   PetscCall(MatDestroy(&A));
50   PetscCall(PetscFinalize());
51   return 0;
52 }
53 
54 /*TEST
55 
56    test:
57 
58    test:
59       suffix: 2
60       args: -baij
61 
62 TEST*/
63