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().\n\n";

 5:  #include petscmat.h

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

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

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

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

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

 48:   VecDuplicate(u,&x);
 49:   VecDuplicate(u,&b);

 51:   MatMultTranspose(C,u,b);

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

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

 61:   VecAXPY(&mone,u,x);
 62:   VecNorm(x,NORM_2,&norm);
 63:   PetscPrintf(PETSC_COMM_SELF,"Norm of error %g\n",norm);

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