Actual source code: ex17.c

petsc-3.7.3 2016-08-01
Report Typos and Errors
  2: static char help[] = "Tests the use of MatSolveTranspose().\n\n";

  4: #include <petscmat.h>

  8: int main(int argc,char **args)
  9: {
 10:   Mat            C,A;
 11:   PetscInt       i,j,m = 5,n = 5,Ii,J;
 13:   PetscScalar    v,five = 5.0,one = 1.0;
 14:   IS             isrow,row,col;
 15:   Vec            x,u,b;
 16:   PetscReal      norm;
 17:   MatFactorInfo  info;

 19:   PetscInitialize(&argc,&args,(char*)0,help);
 20:   PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
 21:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);

 23:   MatCreateSeqAIJ(PETSC_COMM_SELF,m*n,m*n,5,NULL,&C);
 24:   MatSetUp(C);

 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;  Ii = j + n*i;
 30:       if (i>0)   {J = Ii - n; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
 31:       if (i<m-1) {J = Ii + n; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
 32:       if (j>0)   {J = Ii - 1; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
 33:       if (j<n-1) {J = Ii + 1; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
 34:       v = 4.0; MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES);
 35:     }
 36:   }
 37:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 38:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

 40:   ISCreateStride(PETSC_COMM_SELF,(m*n)/2,0,2,&isrow);
 41:   MatZeroRowsIS(C,isrow,five,0,0);

 43:   VecCreateSeq(PETSC_COMM_SELF,m*n,&u);
 44:   VecDuplicate(u,&x);
 45:   VecDuplicate(u,&b);
 46:   VecSet(u,one);

 48:   MatMultTranspose(C,u,b);

 50:   /* Set default ordering to be Quotient Minimum Degree; also read
 51:      orderings from the options database */
 52:   MatGetOrdering(C,MATORDERINGQMD,&row,&col);

 54:   MatFactorInfoInitialize(&info);
 55:   MatGetFactor(C,MATSOLVERPETSC,MAT_FACTOR_LU,&A);
 56:   MatLUFactorSymbolic(A,C,row,col,&info);
 57:   MatLUFactorNumeric(A,C,&info);
 58:   MatSolveTranspose(A,b,x);

 60:   ISView(row,PETSC_VIEWER_STDOUT_SELF);
 61:   VecAXPY(x,-1.0,u);
 62:   VecNorm(x,NORM_2,&norm);
 63:   PetscPrintf(PETSC_COMM_SELF,"Norm of error %g\n",(double)norm);

 65:   ISDestroy(&row);
 66:   ISDestroy(&col);
 67:   ISDestroy(&isrow);
 68:   VecDestroy(&u);
 69:   VecDestroy(&x);
 70:   VecDestroy(&b);
 71:   MatDestroy(&C);
 72:   MatDestroy(&A);
 73:   PetscFinalize();
 74:   return 0;
 75: }