Actual source code: ex117.c

  2: static char help[] = "Tests Cholesky factorization and Matview() for a SBAIJ matrix, (bs=2).\n";
  3: /*
  4:   This code is modified from the code contributed by JUNWANG@uwm.edu on Apr 13, 2007
  5: */

 7:  #include petscmat.h

 11: int main(int argc,char **args)
 12: {
 14:   Mat            mat,fact,B;
 15:   PetscInt       ind1[2],ind2[2];
 16:   PetscScalar    temp[2*2];
 17:   PetscInt       nnz[3];
 18:   IS             perm,colp;
 19:   MatFactorInfo  info;

 21:   PetscInitialize(&argc,&args,0,help);
 22:   nnz[0]=2;nnz[1]=1;nnz[2]=1;

 24:   MatCreateSeqSBAIJ(PETSC_COMM_SELF,2,6,6,0,nnz,&mat);
 25:   ind1[0]=0;ind1[1]=1;
 26:   temp[0]=3;temp[1]=2;temp[2]=0;temp[3]=3;
 27:   MatSetValues(mat,2,ind1,2,ind1,temp,INSERT_VALUES);
 28:   ind2[0]=4;ind2[1]=5;
 29:   temp[0]=1;temp[1]=1;temp[2]=2;temp[3]=1;
 30:   MatSetValues(mat,2,ind1,2,ind2,temp,INSERT_VALUES);
 31:   ind1[0]=2;ind1[1]=3;
 32:   temp[0]=4;temp[1]=1;temp[2]=1;temp[3]=5;
 33:   MatSetValues(mat,2,ind1,2,ind1,temp,INSERT_VALUES);
 34:   ind1[0]=4;ind1[1]=5;
 35:   temp[0]=5;temp[1]=1;temp[2]=1;temp[3]=6;
 36:   MatSetValues(mat,2,ind1,2,ind1,temp,INSERT_VALUES);

 38:   MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);
 39:   MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);

 41:   MatDuplicate(mat,MAT_SHARE_NONZERO_PATTERN,&B);
 42:   ind1[0]=0;ind1[1]=1;
 43:   temp[0]=3;temp[1]=2;temp[2]=0;temp[3]=3;
 44:   MatSetValues(mat,2,ind1,2,ind1,temp,INSERT_VALUES);
 45:   ind2[0]=4;ind2[1]=5;
 46:   temp[0]=1;temp[1]=1;temp[2]=2;temp[3]=1;
 47:   MatSetValues(mat,2,ind1,2,ind2,temp,INSERT_VALUES);
 48:   ind1[0]=2;ind1[1]=3;
 49:   temp[0]=4;temp[1]=1;temp[2]=1;temp[3]=5;
 50:   MatSetValues(mat,2,ind1,2,ind1,temp,INSERT_VALUES);
 51:   ind1[0]=4;ind1[1]=5;
 52:   temp[0]=5;temp[1]=1;temp[2]=1;temp[3]=6;
 53:   MatSetValues(mat,2,ind1,2,ind1,temp,INSERT_VALUES);

 55:   MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);
 56:   MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);

 58:   PetscPrintf(PETSC_COMM_WORLD,"mat: \n");
 59:   MatView(mat,PETSC_VIEWER_STDOUT_SELF);

 61:    // begin cholesky factorization
 62:   MatGetOrdering(mat,MATORDERING_NATURAL,&perm,&colp);
 63:   ISDestroy(colp);
 64:   info.fill=1.0;
 65:   MatGetFactor(mat,MAT_SOLVER_PETSC,MAT_FACTOR_CHOLESKY,&fact);
 66:   MatCholeskyFactorSymbolic(fact,mat,perm,&info);
 67:   MatCholeskyFactorNumeric(fact,mat,&info);
 68:   PetscPrintf(PETSC_COMM_WORLD,"Chol factor: \n");
 69:   MatView(fact, PETSC_VIEWER_STDOUT_SELF);

 71:   ISDestroy(perm);
 72:   MatDestroy(mat);
 73:   MatDestroy(fact);
 74:   MatDestroy(B);
 75:   PetscFinalize();
 76:   return 0;
 77: }