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