Actual source code: ex15.c

  1: /*$Id: ex15.c,v 1.22 2001/08/07 21:30:08 bsmith Exp $*/

  3: static char help[] = "Tests MatNorm(), MatLUFactor(), MatSolve() and MatSolveAdd().nn";

 5:  #include petscmat.h

  7: int main(int argc,char **args)
  8: {
  9:   Mat         C;
 10:   int         i,j,m = 3,n = 3,I,J,ierr;
 11:   PetscTruth  flg;
 12:   PetscScalar v,mone = -1.0,one = 1.0,alpha = 2.0;
 13:   IS          perm,iperm;
 14:   Vec         x,u,b,y;
 15:   PetscReal   norm;

 17:   PetscInitialize(&argc,&args,(char *)0,help);

 19:   MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n,&C);
 20:   MatSetFromOptions(C);
 21:   PetscOptionsHasName(PETSC_NULL,"-symmetric",&flg);
 22:   if (flg) {  /* Treat matrix as symmetric only if we set this flag */
 23:     MatSetOption(C,MAT_SYMMETRIC);
 24:   }

 26:   /* Create the matrix for the five point stencil, YET AGAIN */
 27:   for (i=0; i<m; i++) {
 28:     for (j=0; j<n; j++) {
 29:       v = -1.0;  I = j + n*i;
 30:       if (i>0)   {J = I - n; MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
 31:       if (i<m-1) {J = I + n; MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
 32:       if (j>0)   {J = I - 1; MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
 33:       if (j<n-1) {J = I + 1; MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
 34:       v = 4.0; MatSetValues(C,1,&I,1,&I,&v,INSERT_VALUES);
 35:     }
 36:   }
 37:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 38:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
 39:   MatGetOrdering(C,MATORDERING_RCM,&perm,&iperm);
 40:   MatView(C,PETSC_VIEWER_STDOUT_WORLD);
 41:   ISView(perm,PETSC_VIEWER_STDOUT_SELF);
 42:   VecCreateSeq(PETSC_COMM_SELF,m*n,&u);
 43:   VecSet(&one,u);
 44:   VecDuplicate(u,&x);
 45:   VecDuplicate(u,&b);
 46:   VecDuplicate(u,&y);
 47:   MatMult(C,u,b);
 48:   VecCopy(b,y);
 49:   VecScale(&alpha,y);

 51:   MatNorm(C,NORM_FROBENIUS,&norm);
 52:   PetscPrintf(PETSC_COMM_SELF,"Frobenius norm of matrix %gn",norm);
 53:   MatNorm(C,NORM_1,&norm);
 54:   PetscPrintf(PETSC_COMM_SELF,"One  norm of matrix %gn",norm);
 55:   MatNorm(C,NORM_INFINITY,&norm);
 56:   PetscPrintf(PETSC_COMM_SELF,"Infinity norm of matrix %gn",norm);

 58:   MatLUFactor(C,perm,iperm,PETSC_NULL);
 59:   MatView(C,PETSC_VIEWER_STDOUT_WORLD);

 61:   /* Test MatSolve */
 62:   MatSolve(C,b,x);
 63:   VecView(b,PETSC_VIEWER_STDOUT_SELF);
 64:   VecView(x,PETSC_VIEWER_STDOUT_SELF);
 65:   VecAXPY(&mone,u,x);
 66:   VecNorm(x,NORM_2,&norm);
 67:   PetscPrintf(PETSC_COMM_SELF,"Norm of error %An",norm);

 69:   /* Test MatSolveAdd */
 70:   MatSolveAdd(C,b,y,x);

 72:   VecAXPY(&mone,y,x);
 73:   VecAXPY(&mone,u,x);
 74:   VecNorm(x,NORM_2,&norm);

 76:   PetscPrintf(PETSC_COMM_SELF,"Norm of error %An",norm);

 78:   ISDestroy(perm);
 79:   ISDestroy(iperm);
 80:   VecDestroy(u);
 81:   VecDestroy(y);
 82:   VecDestroy(b);
 83:   VecDestroy(x);
 84:   MatDestroy(C);
 85:   PetscFinalize();
 86:   return 0;
 87: }