xref: /petsc/src/mat/tests/ex15.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests MatNorm(), MatLUFactor(), MatSolve() and MatSolveAdd().\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;
9c4762a1bSJed Brown   PetscInt      i, j, m = 3, n = 3, Ii, J;
10c4762a1bSJed Brown   PetscBool     flg;
11c4762a1bSJed Brown   PetscScalar   v;
12c4762a1bSJed Brown   IS            perm, iperm;
13c4762a1bSJed Brown   Vec           x, u, b, y;
14c4762a1bSJed Brown   PetscReal     norm, tol = PETSC_SMALL;
15c4762a1bSJed Brown   MatFactorInfo info;
16c4762a1bSJed Brown   PetscMPIInt   size;
17c4762a1bSJed Brown 
18327415f7SBarry Smith   PetscFunctionBeginUser;
199566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
209566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
21be096a46SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
229566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &C));
239566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, m * n, m * n));
249566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(C));
259566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
269566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-symmetric", &flg));
27c4762a1bSJed Brown   if (flg) { /* Treat matrix as symmetric only if we set this flag */
289566063dSJacob Faibussowitsch     PetscCall(MatSetOption(C, MAT_SYMMETRIC, PETSC_TRUE));
299566063dSJacob Faibussowitsch     PetscCall(MatSetOption(C, MAT_SYMMETRY_ETERNAL, PETSC_TRUE));
30c4762a1bSJed Brown   }
31c4762a1bSJed Brown 
32c4762a1bSJed Brown   /* Create the matrix for the five point stencil, YET AGAIN */
33c4762a1bSJed Brown   for (i = 0; i < m; i++) {
34c4762a1bSJed Brown     for (j = 0; j < n; j++) {
359371c9d4SSatish Balay       v  = -1.0;
369371c9d4SSatish Balay       Ii = j + n * i;
379371c9d4SSatish Balay       if (i > 0) {
389371c9d4SSatish Balay         J = Ii - n;
399371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
409371c9d4SSatish Balay       }
419371c9d4SSatish Balay       if (i < m - 1) {
429371c9d4SSatish Balay         J = Ii + n;
439371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
449371c9d4SSatish Balay       }
459371c9d4SSatish Balay       if (j > 0) {
469371c9d4SSatish Balay         J = Ii - 1;
479371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
489371c9d4SSatish Balay       }
499371c9d4SSatish Balay       if (j < n - 1) {
509371c9d4SSatish Balay         J = Ii + 1;
519371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
529371c9d4SSatish Balay       }
539371c9d4SSatish Balay       v = 4.0;
549371c9d4SSatish Balay       PetscCall(MatSetValues(C, 1, &Ii, 1, &Ii, &v, INSERT_VALUES));
55c4762a1bSJed Brown     }
56c4762a1bSJed Brown   }
579566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
589566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
599566063dSJacob Faibussowitsch   PetscCall(MatGetOrdering(C, MATORDERINGRCM, &perm, &iperm));
609566063dSJacob Faibussowitsch   PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD));
619566063dSJacob Faibussowitsch   PetscCall(ISView(perm, PETSC_VIEWER_STDOUT_SELF));
629566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, m * n, &u));
639566063dSJacob Faibussowitsch   PetscCall(VecSet(u, 1.0));
649566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &x));
659566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &b));
669566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &y));
679566063dSJacob Faibussowitsch   PetscCall(MatMult(C, u, b));
689566063dSJacob Faibussowitsch   PetscCall(VecCopy(b, y));
699566063dSJacob Faibussowitsch   PetscCall(VecScale(y, 2.0));
70c4762a1bSJed Brown 
719566063dSJacob Faibussowitsch   PetscCall(MatNorm(C, NORM_FROBENIUS, &norm));
729566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "Frobenius norm of matrix %g\n", (double)norm));
739566063dSJacob Faibussowitsch   PetscCall(MatNorm(C, NORM_1, &norm));
749566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "One  norm of matrix %g\n", (double)norm));
759566063dSJacob Faibussowitsch   PetscCall(MatNorm(C, NORM_INFINITY, &norm));
769566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "Infinity norm of matrix %g\n", (double)norm));
77c4762a1bSJed Brown 
789566063dSJacob Faibussowitsch   PetscCall(MatFactorInfoInitialize(&info));
79c4762a1bSJed Brown   info.fill          = 2.0;
80c4762a1bSJed Brown   info.dtcol         = 0.0;
81c4762a1bSJed Brown   info.zeropivot     = 1.e-14;
82c4762a1bSJed Brown   info.pivotinblocks = 1.0;
83c4762a1bSJed Brown 
849566063dSJacob Faibussowitsch   PetscCall(MatLUFactor(C, perm, iperm, &info));
85c4762a1bSJed Brown 
86c4762a1bSJed Brown   /* Test MatSolve */
879566063dSJacob Faibussowitsch   PetscCall(MatSolve(C, b, x));
889566063dSJacob Faibussowitsch   PetscCall(VecView(b, PETSC_VIEWER_STDOUT_SELF));
899566063dSJacob Faibussowitsch   PetscCall(VecView(x, PETSC_VIEWER_STDOUT_SELF));
909566063dSJacob Faibussowitsch   PetscCall(VecAXPY(x, -1.0, u));
919566063dSJacob Faibussowitsch   PetscCall(VecNorm(x, NORM_2, &norm));
9248a46eb9SPierre Jolivet   if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "MatSolve: Norm of error %g\n", (double)norm));
93c4762a1bSJed Brown 
94c4762a1bSJed Brown   /* Test MatSolveAdd */
959566063dSJacob Faibussowitsch   PetscCall(MatSolveAdd(C, b, y, x));
969566063dSJacob Faibussowitsch   PetscCall(VecAXPY(x, -1.0, y));
979566063dSJacob Faibussowitsch   PetscCall(VecAXPY(x, -1.0, u));
989566063dSJacob Faibussowitsch   PetscCall(VecNorm(x, NORM_2, &norm));
9948a46eb9SPierre Jolivet   if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_SELF, "MatSolveAdd(): Norm of error %g\n", (double)norm));
100c4762a1bSJed Brown 
1019566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&perm));
1029566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iperm));
1039566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
1049566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
1059566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b));
1069566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1079566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
1089566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
109b122ec5aSJacob Faibussowitsch   return 0;
110c4762a1bSJed Brown }
111c4762a1bSJed Brown 
112c4762a1bSJed Brown /*TEST
113c4762a1bSJed Brown 
114c4762a1bSJed Brown    test:
115c4762a1bSJed Brown 
116c4762a1bSJed Brown TEST*/
117