1 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 */ 6 7 #include <petscmat.h> 8 9 int main(int argc,char **args) 10 { 11 PetscErrorCode ierr; 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; 19 20 ierr = PetscInitialize(&argc,&args,0,help);if (ierr) return ierr; 21 CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 22 PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!"); 23 24 nnz[0]=2;nnz[1]=1;nnz[2]=1; 25 CHKERRQ(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 CHKERRQ(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 CHKERRQ(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 CHKERRQ(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 CHKERRQ(MatSetValues(mat,2,ind1,2,ind1,temp,INSERT_VALUES)); 39 40 CHKERRQ(MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY)); 41 CHKERRQ(MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY)); 42 43 CHKERRQ(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 CHKERRQ(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 CHKERRQ(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 CHKERRQ(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 CHKERRQ(MatSetValues(mat,2,ind1,2,ind1,temp,INSERT_VALUES)); 56 57 CHKERRQ(MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY)); 58 CHKERRQ(MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY)); 59 60 CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"mat: \n")); 61 CHKERRQ(MatView(mat,PETSC_VIEWER_STDOUT_SELF)); 62 63 /* begin cholesky factorization */ 64 CHKERRQ(MatGetOrdering(mat,MATORDERINGNATURAL,&perm,&colp)); 65 CHKERRQ(ISDestroy(&colp)); 66 67 info.fill=1.0; 68 CHKERRQ(MatGetFactor(mat,MATSOLVERPETSC,MAT_FACTOR_CHOLESKY,&fact)); 69 CHKERRQ(MatCholeskyFactorSymbolic(fact,mat,perm,&info)); 70 CHKERRQ(MatCholeskyFactorNumeric(fact,mat,&info)); 71 72 CHKERRQ(ISDestroy(&perm)); 73 CHKERRQ(MatDestroy(&mat)); 74 CHKERRQ(MatDestroy(&fact)); 75 CHKERRQ(MatDestroy(&B)); 76 ierr = PetscFinalize(); 77 return ierr; 78 } 79