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