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