Actual source code: ex15.c


  2: static char help[] = "Tests MatNorm(), MatLUFactor(), MatSolve() and MatSolveAdd().\n\n";

  4: #include <petscmat.h>

  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;

 18:   PetscInitialize(&argc,&args,(char*)0,help);
 19:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 21:   MatCreate(PETSC_COMM_WORLD,&C);
 22:   MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);
 23:   MatSetFromOptions(C);
 24:   MatSetUp(C);
 25:   PetscOptionsHasName(NULL,NULL,"-symmetric",&flg);
 26:   if (flg) {  /* Treat matrix as symmetric only if we set this flag */
 27:     MatSetOption(C,MAT_SYMMETRIC,PETSC_TRUE);
 28:     MatSetOption(C,MAT_SYMMETRY_ETERNAL,PETSC_TRUE);
 29:   }

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

 56:   MatNorm(C,NORM_FROBENIUS,&norm);
 57:   PetscPrintf(PETSC_COMM_SELF,"Frobenius norm of matrix %g\n",(double)norm);
 58:   MatNorm(C,NORM_1,&norm);
 59:   PetscPrintf(PETSC_COMM_SELF,"One  norm of matrix %g\n",(double)norm);
 60:   MatNorm(C,NORM_INFINITY,&norm);
 61:   PetscPrintf(PETSC_COMM_SELF,"Infinity norm of matrix %g\n",(double)norm);

 63:   MatFactorInfoInitialize(&info);
 64:   info.fill          = 2.0;
 65:   info.dtcol         = 0.0;
 66:   info.zeropivot     = 1.e-14;
 67:   info.pivotinblocks = 1.0;

 69:   MatLUFactor(C,perm,iperm,&info);

 71:   /* Test MatSolve */
 72:   MatSolve(C,b,x);
 73:   VecView(b,PETSC_VIEWER_STDOUT_SELF);
 74:   VecView(x,PETSC_VIEWER_STDOUT_SELF);
 75:   VecAXPY(x,-1.0,u);
 76:   VecNorm(x,NORM_2,&norm);
 77:   if (norm > tol) {
 78:     PetscPrintf(PETSC_COMM_SELF,"MatSolve: Norm of error %g\n",(double)norm);
 79:   }

 81:   /* Test MatSolveAdd */
 82:   MatSolveAdd(C,b,y,x);
 83:   VecAXPY(x,-1.0,y);
 84:   VecAXPY(x,-1.0,u);
 85:   VecNorm(x,NORM_2,&norm);
 86:   if (norm > tol) {
 87:     PetscPrintf(PETSC_COMM_SELF,"MatSolveAdd(): Norm of error %g\n",(double)norm);
 88:   }

 90:   ISDestroy(&perm);
 91:   ISDestroy(&iperm);
 92:   VecDestroy(&u);
 93:   VecDestroy(&y);
 94:   VecDestroy(&b);
 95:   VecDestroy(&x);
 96:   MatDestroy(&C);
 97:   PetscFinalize();
 98:   return 0;
 99: }

101: /*TEST

103:    test:

105: TEST*/