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 ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 22 if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!"); 23 24 nnz[0]=2;nnz[1]=1;nnz[2]=1; 25 ierr = MatCreateSeqSBAIJ(PETSC_COMM_SELF,2,6,6,0,nnz,&mat);CHKERRQ(ierr); 26 27 ind1[0]=0;ind1[1]=1; 28 temp[0]=3;temp[1]=2;temp[2]=0;temp[3]=3; 29 ierr = MatSetValues(mat,2,ind1,2,ind1,temp,INSERT_VALUES);CHKERRQ(ierr); 30 ind2[0]=4;ind2[1]=5; 31 temp[0]=1;temp[1]=1;temp[2]=2;temp[3]=1; 32 ierr = MatSetValues(mat,2,ind1,2,ind2,temp,INSERT_VALUES);CHKERRQ(ierr); 33 ind1[0]=2;ind1[1]=3; 34 temp[0]=4;temp[1]=1;temp[2]=1;temp[3]=5; 35 ierr = MatSetValues(mat,2,ind1,2,ind1,temp,INSERT_VALUES);CHKERRQ(ierr); 36 ind1[0]=4;ind1[1]=5; 37 temp[0]=5;temp[1]=1;temp[2]=1;temp[3]=6; 38 ierr = MatSetValues(mat,2,ind1,2,ind1,temp,INSERT_VALUES);CHKERRQ(ierr); 39 40 ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 41 ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 42 43 ierr = MatDuplicate(mat,MAT_SHARE_NONZERO_PATTERN,&B);CHKERRQ(ierr); 44 ind1[0]=0;ind1[1]=1; 45 temp[0]=3;temp[1]=2;temp[2]=0;temp[3]=3; 46 ierr = MatSetValues(mat,2,ind1,2,ind1,temp,INSERT_VALUES);CHKERRQ(ierr); 47 ind2[0]=4;ind2[1]=5; 48 temp[0]=1;temp[1]=1;temp[2]=2;temp[3]=1; 49 ierr = MatSetValues(mat,2,ind1,2,ind2,temp,INSERT_VALUES);CHKERRQ(ierr); 50 ind1[0]=2;ind1[1]=3; 51 temp[0]=4;temp[1]=1;temp[2]=1;temp[3]=5; 52 ierr = MatSetValues(mat,2,ind1,2,ind1,temp,INSERT_VALUES);CHKERRQ(ierr); 53 ind1[0]=4;ind1[1]=5; 54 temp[0]=5;temp[1]=1;temp[2]=1;temp[3]=6; 55 ierr = MatSetValues(mat,2,ind1,2,ind1,temp,INSERT_VALUES);CHKERRQ(ierr); 56 57 ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 58 ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 59 60 ierr = PetscPrintf(PETSC_COMM_WORLD,"mat: \n");CHKERRQ(ierr); 61 ierr = MatView(mat,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 62 63 /* begin cholesky factorization */ 64 ierr = MatGetOrdering(mat,MATORDERINGNATURAL,&perm,&colp);CHKERRQ(ierr); 65 ierr = ISDestroy(&colp);CHKERRQ(ierr); 66 67 info.fill=1.0; 68 ierr = MatGetFactor(mat,MATSOLVERPETSC,MAT_FACTOR_CHOLESKY,&fact);CHKERRQ(ierr); 69 ierr = MatCholeskyFactorSymbolic(fact,mat,perm,&info);CHKERRQ(ierr); 70 ierr = MatCholeskyFactorNumeric(fact,mat,&info);CHKERRQ(ierr); 71 72 ierr = ISDestroy(&perm);CHKERRQ(ierr); 73 ierr = MatDestroy(&mat);CHKERRQ(ierr); 74 ierr = MatDestroy(&fact);CHKERRQ(ierr); 75 ierr = MatDestroy(&B);CHKERRQ(ierr); 76 ierr = PetscFinalize(); 77 return ierr; 78 } 79