xref: /petsc/src/mat/tests/ex13.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests copying and ordering uniprocessor row-based sparse matrices.\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6*d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
7*d71ae5a4SJacob Faibussowitsch {
8c4762a1bSJed Brown   Mat         C, A;
9c4762a1bSJed Brown   PetscInt    i, j, m = 5, n = 5, Ii, J;
10c4762a1bSJed Brown   PetscScalar v;
11c4762a1bSJed Brown   IS          perm, iperm;
12c4762a1bSJed Brown 
13327415f7SBarry Smith   PetscFunctionBeginUser;
149566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
159566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, m * n, m * n, 5, NULL, &C));
16c4762a1bSJed Brown   /* create the matrix for the five point stencil, YET AGAIN*/
17c4762a1bSJed Brown   for (i = 0; i < m; i++) {
18c4762a1bSJed Brown     for (j = 0; j < n; j++) {
199371c9d4SSatish Balay       v  = -1.0;
209371c9d4SSatish Balay       Ii = j + n * i;
219371c9d4SSatish Balay       if (i > 0) {
229371c9d4SSatish Balay         J = Ii - n;
239371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
249371c9d4SSatish Balay       }
259371c9d4SSatish Balay       if (i < m - 1) {
269371c9d4SSatish Balay         J = Ii + n;
279371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
289371c9d4SSatish Balay       }
299371c9d4SSatish Balay       if (j > 0) {
309371c9d4SSatish Balay         J = Ii - 1;
319371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
329371c9d4SSatish Balay       }
339371c9d4SSatish Balay       if (j < n - 1) {
349371c9d4SSatish Balay         J = Ii + 1;
359371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
369371c9d4SSatish Balay       }
379371c9d4SSatish Balay       v = 4.0;
389371c9d4SSatish Balay       PetscCall(MatSetValues(C, 1, &Ii, 1, &Ii, &v, INSERT_VALUES));
39c4762a1bSJed Brown     }
40c4762a1bSJed Brown   }
419566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
429566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
43c4762a1bSJed Brown 
449566063dSJacob Faibussowitsch   PetscCall(MatConvert(C, MATSAME, MAT_INITIAL_MATRIX, &A));
45c4762a1bSJed Brown 
469566063dSJacob Faibussowitsch   PetscCall(MatGetOrdering(A, MATORDERINGND, &perm, &iperm));
479566063dSJacob Faibussowitsch   PetscCall(ISView(perm, PETSC_VIEWER_STDOUT_SELF));
489566063dSJacob Faibussowitsch   PetscCall(ISView(iperm, PETSC_VIEWER_STDOUT_SELF));
499566063dSJacob Faibussowitsch   PetscCall(MatView(A, PETSC_VIEWER_STDOUT_SELF));
50c4762a1bSJed Brown 
519566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&perm));
529566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iperm));
539566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
549566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
559566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
56b122ec5aSJacob Faibussowitsch   return 0;
57c4762a1bSJed Brown }
58c4762a1bSJed Brown 
59c4762a1bSJed Brown /*TEST
60c4762a1bSJed Brown 
61c4762a1bSJed Brown    test:
628cc725e6SPierre Jolivet       filter: grep -v " MPI process"
63c4762a1bSJed Brown 
64c4762a1bSJed Brown TEST*/
65