xref: /petsc/src/mat/tests/ex61.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests MatSeq(B)AIJSetColumnIndices().\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown /*
7c4762a1bSJed Brown       Generate the following matrix:
8c4762a1bSJed Brown 
9c4762a1bSJed Brown          1 0 3
10c4762a1bSJed Brown          1 2 3
11c4762a1bSJed Brown          0 0 3
12c4762a1bSJed Brown */
13*d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
14*d71ae5a4SJacob Faibussowitsch {
15c4762a1bSJed Brown   Mat         A;
16c4762a1bSJed Brown   PetscScalar v;
17c4762a1bSJed Brown   PetscInt    i, j, rowlens[] = {2, 3, 1}, cols[] = {0, 2, 0, 1, 2, 2};
18c4762a1bSJed Brown   PetscBool   flg;
19c4762a1bSJed Brown 
20327415f7SBarry Smith   PetscFunctionBeginUser;
219566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-baij", &flg));
23c4762a1bSJed Brown   if (flg) {
249566063dSJacob Faibussowitsch     PetscCall(MatCreateSeqBAIJ(PETSC_COMM_WORLD, 1, 3, 3, 0, rowlens, &A));
259566063dSJacob Faibussowitsch     PetscCall(MatSeqBAIJSetColumnIndices(A, cols));
26c4762a1bSJed Brown   } else {
279566063dSJacob Faibussowitsch     PetscCall(MatCreateSeqAIJ(PETSC_COMM_WORLD, 3, 3, 0, rowlens, &A));
289566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetColumnIndices(A, cols));
29c4762a1bSJed Brown   }
30c4762a1bSJed Brown 
319371c9d4SSatish Balay   i = 0;
329371c9d4SSatish Balay   j = 0;
339371c9d4SSatish Balay   v = 1.0;
349566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &i, 1, &j, &v, INSERT_VALUES));
359371c9d4SSatish Balay   i = 0;
369371c9d4SSatish Balay   j = 2;
379371c9d4SSatish Balay   v = 3.0;
389566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &i, 1, &j, &v, INSERT_VALUES));
39c4762a1bSJed Brown 
409371c9d4SSatish Balay   i = 1;
419371c9d4SSatish Balay   j = 0;
429371c9d4SSatish Balay   v = 1.0;
439566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &i, 1, &j, &v, INSERT_VALUES));
449371c9d4SSatish Balay   i = 1;
459371c9d4SSatish Balay   j = 1;
469371c9d4SSatish Balay   v = 2.0;
479566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &i, 1, &j, &v, INSERT_VALUES));
489371c9d4SSatish Balay   i = 1;
499371c9d4SSatish Balay   j = 2;
509371c9d4SSatish Balay   v = 3.0;
519566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &i, 1, &j, &v, INSERT_VALUES));
52c4762a1bSJed Brown 
539371c9d4SSatish Balay   i = 2;
549371c9d4SSatish Balay   j = 2;
559371c9d4SSatish Balay   v = 3.0;
569566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &i, 1, &j, &v, INSERT_VALUES));
57c4762a1bSJed Brown 
589566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
599566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
609566063dSJacob Faibussowitsch   PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD));
61c4762a1bSJed Brown 
629566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
639566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
64b122ec5aSJacob Faibussowitsch   return 0;
65c4762a1bSJed Brown }
66c4762a1bSJed Brown 
67c4762a1bSJed Brown /*TEST
68c4762a1bSJed Brown 
69c4762a1bSJed Brown    test:
70c4762a1bSJed Brown 
71c4762a1bSJed Brown    test:
72c4762a1bSJed Brown       suffix: 2
73c4762a1bSJed Brown       args: -baij
74c4762a1bSJed Brown 
75c4762a1bSJed Brown TEST*/
76