xref: /petsc/src/mat/tests/ex61.c (revision 732aec7a18f2199fb53bb9a2f3aef439a834ce31) !
1c4762a1bSJed Brown static char help[] = "Tests MatSeq(B)AIJSetColumnIndices().\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown 
5c4762a1bSJed Brown /*
6c4762a1bSJed Brown       Generate the following matrix:
7c4762a1bSJed Brown 
8c4762a1bSJed Brown          1 0 3
9c4762a1bSJed Brown          1 2 3
10c4762a1bSJed Brown          0 0 3
11c4762a1bSJed Brown */
main(int argc,char ** args)12d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
13d71ae5a4SJacob Faibussowitsch {
14c4762a1bSJed Brown   Mat         A;
15c4762a1bSJed Brown   PetscScalar v;
16c4762a1bSJed Brown   PetscInt    i, j, rowlens[] = {2, 3, 1}, cols[] = {0, 2, 0, 1, 2, 2};
17c4762a1bSJed Brown   PetscBool   flg;
18c4762a1bSJed Brown 
19327415f7SBarry Smith   PetscFunctionBeginUser;
20*c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &args, NULL, help));
219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-baij", &flg));
22c4762a1bSJed Brown   if (flg) {
239566063dSJacob Faibussowitsch     PetscCall(MatCreateSeqBAIJ(PETSC_COMM_WORLD, 1, 3, 3, 0, rowlens, &A));
249566063dSJacob Faibussowitsch     PetscCall(MatSeqBAIJSetColumnIndices(A, cols));
25c4762a1bSJed Brown   } else {
269566063dSJacob Faibussowitsch     PetscCall(MatCreateSeqAIJ(PETSC_COMM_WORLD, 3, 3, 0, rowlens, &A));
279566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetColumnIndices(A, cols));
28c4762a1bSJed Brown   }
29c4762a1bSJed Brown 
309371c9d4SSatish Balay   i = 0;
319371c9d4SSatish Balay   j = 0;
329371c9d4SSatish Balay   v = 1.0;
339566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &i, 1, &j, &v, INSERT_VALUES));
349371c9d4SSatish Balay   i = 0;
359371c9d4SSatish Balay   j = 2;
369371c9d4SSatish Balay   v = 3.0;
379566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &i, 1, &j, &v, INSERT_VALUES));
38c4762a1bSJed Brown 
399371c9d4SSatish Balay   i = 1;
409371c9d4SSatish Balay   j = 0;
419371c9d4SSatish Balay   v = 1.0;
429566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &i, 1, &j, &v, INSERT_VALUES));
439371c9d4SSatish Balay   i = 1;
449371c9d4SSatish Balay   j = 1;
459371c9d4SSatish Balay   v = 2.0;
469566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &i, 1, &j, &v, INSERT_VALUES));
479371c9d4SSatish Balay   i = 1;
489371c9d4SSatish Balay   j = 2;
499371c9d4SSatish Balay   v = 3.0;
509566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &i, 1, &j, &v, INSERT_VALUES));
51c4762a1bSJed Brown 
529371c9d4SSatish Balay   i = 2;
539371c9d4SSatish Balay   j = 2;
549371c9d4SSatish Balay   v = 3.0;
559566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &i, 1, &j, &v, INSERT_VALUES));
56c4762a1bSJed Brown 
579566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
589566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
599566063dSJacob Faibussowitsch   PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD));
60c4762a1bSJed Brown 
619566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
629566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
63b122ec5aSJacob Faibussowitsch   return 0;
64c4762a1bSJed Brown }
65c4762a1bSJed Brown 
66c4762a1bSJed Brown /*TEST
67c4762a1bSJed Brown 
68c4762a1bSJed Brown    test:
69c4762a1bSJed Brown 
70c4762a1bSJed Brown    test:
71c4762a1bSJed Brown       suffix: 2
72c4762a1bSJed Brown       args: -baij
73c4762a1bSJed Brown 
74c4762a1bSJed Brown TEST*/
75