Actual source code: ex117.c

petsc-3.13.6 2020-09-29
Report Typos and Errors

  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: {
 12:   Mat            mat,fact,B;
 13:   PetscInt       ind1[2],ind2[2];
 14:   PetscScalar    temp[4];
 15:   PetscInt       nnz[3];
 16:   IS             perm,colp;
 17:   MatFactorInfo  info;
 18:   PetscMPIInt    size;

 20:   PetscInitialize(&argc,&args,0,help);if (ierr) return ierr;
 21:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 22:   if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!");

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

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

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

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

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

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

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

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