Actual source code: ex17.c

petsc-3.13.6 2020-09-29
Report Typos and Errors

  2: static char help[] = "Tests the use of MatSolveTranspose().\n\n";

  4:  #include <petscmat.h>

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

 17:   PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
 18:   PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
 19:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);

 21:   MatCreateSeqAIJ(PETSC_COMM_SELF,m*n,m*n,5,NULL,&C);
 22:   MatSetUp(C);

 24:   /* create the matrix for the five point stencil, YET AGAIN*/
 25:   for (i=0; i<m; i++) {
 26:     for (j=0; j<n; j++) {
 27:       v = -1.0;  Ii = j + n*i;
 28:       if (i>0)   {J = Ii - n; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
 29:       if (i<m-1) {J = Ii + n; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
 30:       if (j>0)   {J = Ii - 1; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
 31:       if (j<n-1) {J = Ii + 1; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
 32:       v = 4.0; MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES);
 33:     }
 34:   }
 35:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 36:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

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

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

 46:   MatMultTranspose(C,u,b);

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

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

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

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


 76: /*TEST

 78:    test:

 80: TEST*/