xref: /petsc/src/mat/tests/ex15.c (revision ebead697dbf761eb322f829370bbe90b3bd93fa3)
1 
2 static char help[] = "Tests MatNorm(), MatLUFactor(), MatSolve() and MatSolveAdd().\n\n";
3 
4 #include <petscmat.h>
5 
6 int main(int argc,char **args)
7 {
8   Mat            C;
9   PetscInt       i,j,m = 3,n = 3,Ii,J;
10   PetscBool      flg;
11   PetscScalar    v;
12   IS             perm,iperm;
13   Vec            x,u,b,y;
14   PetscReal      norm,tol=PETSC_SMALL;
15   MatFactorInfo  info;
16   PetscMPIInt    size;
17 
18   PetscFunctionBeginUser;
19   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
20   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
21   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
22   PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
23   PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n));
24   PetscCall(MatSetFromOptions(C));
25   PetscCall(MatSetUp(C));
26   PetscCall(PetscOptionsHasName(NULL,NULL,"-symmetric",&flg));
27   if (flg) {  /* Treat matrix as symmetric only if we set this flag */
28     PetscCall(MatSetOption(C,MAT_SYMMETRIC,PETSC_TRUE));
29     PetscCall(MatSetOption(C,MAT_SYMMETRY_ETERNAL,PETSC_TRUE));
30   }
31 
32   /* Create the matrix for the five point stencil, YET AGAIN */
33   for (i=0; i<m; i++) {
34     for (j=0; j<n; j++) {
35       v = -1.0;  Ii = j + n*i;
36       if (i>0)   {J = Ii - n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
37       if (i<m-1) {J = Ii + n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
38       if (j>0)   {J = Ii - 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
39       if (j<n-1) {J = Ii + 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
40       v = 4.0; PetscCall(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES));
41     }
42   }
43   PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
44   PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
45   PetscCall(MatGetOrdering(C,MATORDERINGRCM,&perm,&iperm));
46   PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
47   PetscCall(ISView(perm,PETSC_VIEWER_STDOUT_SELF));
48   PetscCall(VecCreateSeq(PETSC_COMM_SELF,m*n,&u));
49   PetscCall(VecSet(u,1.0));
50   PetscCall(VecDuplicate(u,&x));
51   PetscCall(VecDuplicate(u,&b));
52   PetscCall(VecDuplicate(u,&y));
53   PetscCall(MatMult(C,u,b));
54   PetscCall(VecCopy(b,y));
55   PetscCall(VecScale(y,2.0));
56 
57   PetscCall(MatNorm(C,NORM_FROBENIUS,&norm));
58   PetscCall(PetscPrintf(PETSC_COMM_SELF,"Frobenius norm of matrix %g\n",(double)norm));
59   PetscCall(MatNorm(C,NORM_1,&norm));
60   PetscCall(PetscPrintf(PETSC_COMM_SELF,"One  norm of matrix %g\n",(double)norm));
61   PetscCall(MatNorm(C,NORM_INFINITY,&norm));
62   PetscCall(PetscPrintf(PETSC_COMM_SELF,"Infinity norm of matrix %g\n",(double)norm));
63 
64   PetscCall(MatFactorInfoInitialize(&info));
65   info.fill          = 2.0;
66   info.dtcol         = 0.0;
67   info.zeropivot     = 1.e-14;
68   info.pivotinblocks = 1.0;
69 
70   PetscCall(MatLUFactor(C,perm,iperm,&info));
71 
72   /* Test MatSolve */
73   PetscCall(MatSolve(C,b,x));
74   PetscCall(VecView(b,PETSC_VIEWER_STDOUT_SELF));
75   PetscCall(VecView(x,PETSC_VIEWER_STDOUT_SELF));
76   PetscCall(VecAXPY(x,-1.0,u));
77   PetscCall(VecNorm(x,NORM_2,&norm));
78   if (norm > tol) {
79     PetscCall(PetscPrintf(PETSC_COMM_SELF,"MatSolve: Norm of error %g\n",(double)norm));
80   }
81 
82   /* Test MatSolveAdd */
83   PetscCall(MatSolveAdd(C,b,y,x));
84   PetscCall(VecAXPY(x,-1.0,y));
85   PetscCall(VecAXPY(x,-1.0,u));
86   PetscCall(VecNorm(x,NORM_2,&norm));
87   if (norm > tol) {
88     PetscCall(PetscPrintf(PETSC_COMM_SELF,"MatSolveAdd(): Norm of error %g\n",(double)norm));
89   }
90 
91   PetscCall(ISDestroy(&perm));
92   PetscCall(ISDestroy(&iperm));
93   PetscCall(VecDestroy(&u));
94   PetscCall(VecDestroy(&y));
95   PetscCall(VecDestroy(&b));
96   PetscCall(VecDestroy(&x));
97   PetscCall(MatDestroy(&C));
98   PetscCall(PetscFinalize());
99   return 0;
100 }
101 
102 /*TEST
103 
104    test:
105 
106 TEST*/
107