Actual source code: ex62.c

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

  3: static char help[] = "Tests the use of MatSolveTranspose().nn";

 5:  #include petscmat.h

  7: int main(int argc,char **args)
  8: {
  9:   Mat         C,A;
 10:   int         i,j,m,ierr,size;
 11:   IS          row,col;
 12:   Vec         x,u,b;
 13:   PetscReal   norm;
 14:   PetscViewer fd;
 15:   char        type[256];
 16:   char        file[128];
 17:   PetscScalar one = 1.0,mone = -1.0;
 18:   PetscTruth  flg;

 20:   PetscInitialize(&argc,&args,(char *)0,help);
 21:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 22:   if (size > 1) SETERRQ(1,"Can only run on one processor");

 24:   PetscOptionsGetString(PETSC_NULL,"-f",file,127,&flg);
 25:   if (!flg) SETERRQ(1,"Must indicate binary file with the -f option");
 26:   /* 
 27:      Open binary file.  Note that we use PETSC_BINARY_RDONLY to indicate
 28:      reading from this file.
 29:   */
 30:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,PETSC_BINARY_RDONLY,&fd);

 32:   /* 
 33:      Determine matrix format to be used (specified at runtime).
 34:      See the manpage for MatLoad() for available formats.
 35:   */
 36:   PetscStrcpy(type,MATSEQAIJ);
 37:   PetscOptionsGetString(PETSC_NULL,"-mat_type",type,256,PETSC_NULL);

 39:   /*
 40:      Load the matrix and vector; then destroy the viewer.
 41:   */
 42:   MatLoad(fd,type,&C);
 43:   VecLoad(fd,&u);
 44:   PetscViewerDestroy(fd);

 46:   VecDuplicate(u,&x);
 47:   VecDuplicate(u,&b);

 49:   MatMultTranspose(C,u,b);

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

 55:   MatLUFactorSymbolic(C,row,col,PETSC_NULL,&A);
 56:   MatLUFactorNumeric(C,&A);
 57:   MatSolveTranspose(A,b,x);

 59:   VecAXPY(&mone,u,x);
 60:   VecNorm(x,NORM_2,&norm);
 61:   PetscPrintf(PETSC_COMM_SELF,"Norm of error %gn",norm);

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