Actual source code: ex117.c


  2: static char help[] = "Tests Cholesky factorization 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>

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

 19:   PetscInitialize(&argc,&args,0,help);
 20:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

 23:   nnz[0]=2;nnz[1]=1;nnz[2]=1;
 24:   MatCreateSeqSBAIJ(PETSC_COMM_SELF,2,6,6,0,nnz,&mat);

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

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

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

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

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

 62:   /* begin cholesky factorization */
 63:   MatGetOrdering(mat,MATORDERINGNATURAL,&perm,&colp);
 64:   ISDestroy(&colp);

 66:   info.fill=1.0;
 67:   MatGetFactor(mat,MATSOLVERPETSC,MAT_FACTOR_CHOLESKY,&fact);
 68:   MatCholeskyFactorSymbolic(fact,mat,perm,&info);
 69:   MatCholeskyFactorNumeric(fact,mat,&info);

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