Actual source code: ex1.c

  1: /*$Id: ex1.c,v 1.20 2001/03/23 23:22:29 balay Exp $*/

  3: static char help[] = "Tests LU and Cholesky factorization for a dense matrix.nn";

 5:  #include petscmat.h

  7: int main(int argc,char **argv)
  8: {
  9:   Mat        mat,fact;
 10:   MatInfo    info;
 11:   int        m = 10,n = 10,i = 4,ierr,rstart,rend;
 12:   Scalar     value = 1.0;
 13:   Vec        x,y,b;
 14:   double     norm;

 16:   PetscInitialize(&argc,&argv,(char*) 0,help);

 18:   VecCreate(PETSC_COMM_WORLD,PETSC_DECIDE,m,&y);
 19:   VecSetFromOptions(y);
 20:   VecDuplicate(y,&x);
 21:   VecSet(&value,x);
 22:   VecCreate(PETSC_COMM_WORLD,PETSC_DECIDE,n,&b);
 23:   VecSetFromOptions(b);

 25:   MatCreateSeqDense(PETSC_COMM_WORLD,m,n,PETSC_NULL,&mat);

 27:   MatGetOwnershipRange(mat,&rstart,&rend);
 28:   for (i=rstart; i<rend; i++) {
 29:     value = (double)i+1;
 30:     MatSetValues(mat,1,&i,1,&i,&value,INSERT_VALUES);
 31:   }
 32:   MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);
 33:   MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);

 35:   MatGetInfo(mat,MAT_LOCAL,&info);
 36:   PetscPrintf(PETSC_COMM_WORLD,"matrix nonzeros = %d, allocated nonzeros = %dn",
 37:     (int)info.nz_used,(int)info.nz_allocated);

 39:   /* Cholesky factorization is not yet in place for this matrix format */
 40:   MatMult(mat,x,b);
 41:   MatConvert(mat,MATSAME,&fact);
 42:   MatCholeskyFactor(fact,0,1.0);
 43:   MatSolve(fact,b,y);
 44:   MatDestroy(fact);
 45:   value = -1.0; VecAXPY(&value,x,y);
 46:   VecNorm(y,NORM_2,&norm);
 47:   PetscPrintf(PETSC_COMM_WORLD,"Norm of error for Cholesky %An",norm);
 48:   MatCholeskyFactorSymbolic(mat,0,1.0,&fact);
 49:   MatCholeskyFactorNumeric(mat,&fact);
 50:   MatSolve(fact,b,y);
 51:   value = -1.0; VecAXPY(&value,x,y);
 52:   VecNorm(y,NORM_2,&norm);
 53:   PetscPrintf(PETSC_COMM_WORLD,"Norm of error for Cholesky %An",norm);
 54:   MatDestroy(fact);

 56:   i = m-1; value = 1.0;
 57:   MatSetValues(mat,1,&i,1,&i,&value,INSERT_VALUES);
 58:   MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);
 59:   MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);
 60:   MatMult(mat,x,b);
 61:   MatConvert(mat,MATSAME,&fact);
 62:   MatLUFactor(fact,0,0,PETSC_NULL);
 63:   MatSolve(fact,b,y);
 64:   value = -1.0; VecAXPY(&value,x,y);
 65:   VecNorm(y,NORM_2,&norm);
 66:   PetscPrintf(PETSC_COMM_WORLD,"Norm of error for LU %An",norm);
 67:   MatDestroy(fact);

 69:   MatLUFactorSymbolic(mat,0,0,PETSC_NULL,&fact);
 70:   MatLUFactorNumeric(mat,&fact);
 71:   MatSolve(fact,b,y);
 72:   value = -1.0; VecAXPY(&value,x,y);
 73:   VecNorm(y,NORM_2,&norm);
 74:   PetscPrintf(PETSC_COMM_WORLD,"Norm of error for LU %An",norm);
 75:   MatDestroy(fact);
 76:   MatDestroy(mat);
 77:   VecDestroy(x);
 78:   VecDestroy(b);
 79:   VecDestroy(y);

 81:   PetscFinalize();
 82:   return 0;
 83: }
 84: