1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests Cholesky factorization for a SBAIJ matrix, (bs=2).\n"; 3c4762a1bSJed Brown /* 4c4762a1bSJed Brown This code is modified from the code contributed by JUNWANG@uwm.edu on Apr 13, 2007 5c4762a1bSJed Brown */ 6c4762a1bSJed Brown 7c4762a1bSJed Brown #include <petscmat.h> 8c4762a1bSJed Brown 9c4762a1bSJed Brown int main(int argc,char **args) 10c4762a1bSJed Brown { 11c4762a1bSJed Brown Mat mat,fact,B; 12c4762a1bSJed Brown PetscInt ind1[2],ind2[2]; 13c4762a1bSJed Brown PetscScalar temp[4]; 14c4762a1bSJed Brown PetscInt nnz[3]; 15c4762a1bSJed Brown IS perm,colp; 16c4762a1bSJed Brown MatFactorInfo info; 17c4762a1bSJed Brown PetscMPIInt size; 18c4762a1bSJed Brown 19*327415f7SBarry Smith PetscFunctionBeginUser; 209566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,0,help)); 219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 22be096a46SBarry Smith PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 23c4762a1bSJed Brown 24c4762a1bSJed Brown nnz[0]=2;nnz[1]=1;nnz[2]=1; 259566063dSJacob Faibussowitsch PetscCall(MatCreateSeqSBAIJ(PETSC_COMM_SELF,2,6,6,0,nnz,&mat)); 26c4762a1bSJed Brown 27c4762a1bSJed Brown ind1[0]=0;ind1[1]=1; 28c4762a1bSJed Brown temp[0]=3;temp[1]=2;temp[2]=0;temp[3]=3; 299566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat,2,ind1,2,ind1,temp,INSERT_VALUES)); 30c4762a1bSJed Brown ind2[0]=4;ind2[1]=5; 31c4762a1bSJed Brown temp[0]=1;temp[1]=1;temp[2]=2;temp[3]=1; 329566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat,2,ind1,2,ind2,temp,INSERT_VALUES)); 33c4762a1bSJed Brown ind1[0]=2;ind1[1]=3; 34c4762a1bSJed Brown temp[0]=4;temp[1]=1;temp[2]=1;temp[3]=5; 359566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat,2,ind1,2,ind1,temp,INSERT_VALUES)); 36c4762a1bSJed Brown ind1[0]=4;ind1[1]=5; 37c4762a1bSJed Brown temp[0]=5;temp[1]=1;temp[2]=1;temp[3]=6; 389566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat,2,ind1,2,ind1,temp,INSERT_VALUES)); 39c4762a1bSJed Brown 409566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY)); 419566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY)); 42c4762a1bSJed Brown 439566063dSJacob Faibussowitsch PetscCall(MatDuplicate(mat,MAT_SHARE_NONZERO_PATTERN,&B)); 44c4762a1bSJed Brown ind1[0]=0;ind1[1]=1; 45c4762a1bSJed Brown temp[0]=3;temp[1]=2;temp[2]=0;temp[3]=3; 469566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat,2,ind1,2,ind1,temp,INSERT_VALUES)); 47c4762a1bSJed Brown ind2[0]=4;ind2[1]=5; 48c4762a1bSJed Brown temp[0]=1;temp[1]=1;temp[2]=2;temp[3]=1; 499566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat,2,ind1,2,ind2,temp,INSERT_VALUES)); 50c4762a1bSJed Brown ind1[0]=2;ind1[1]=3; 51c4762a1bSJed Brown temp[0]=4;temp[1]=1;temp[2]=1;temp[3]=5; 529566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat,2,ind1,2,ind1,temp,INSERT_VALUES)); 53c4762a1bSJed Brown ind1[0]=4;ind1[1]=5; 54c4762a1bSJed Brown temp[0]=5;temp[1]=1;temp[2]=1;temp[3]=6; 559566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat,2,ind1,2,ind1,temp,INSERT_VALUES)); 56c4762a1bSJed Brown 579566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY)); 589566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY)); 59c4762a1bSJed Brown 609566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"mat: \n")); 619566063dSJacob Faibussowitsch PetscCall(MatView(mat,PETSC_VIEWER_STDOUT_SELF)); 62c4762a1bSJed Brown 63c4762a1bSJed Brown /* begin cholesky factorization */ 649566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(mat,MATORDERINGNATURAL,&perm,&colp)); 659566063dSJacob Faibussowitsch PetscCall(ISDestroy(&colp)); 66c4762a1bSJed Brown 67c4762a1bSJed Brown info.fill=1.0; 689566063dSJacob Faibussowitsch PetscCall(MatGetFactor(mat,MATSOLVERPETSC,MAT_FACTOR_CHOLESKY,&fact)); 699566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorSymbolic(fact,mat,perm,&info)); 709566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorNumeric(fact,mat,&info)); 71c4762a1bSJed Brown 729566063dSJacob Faibussowitsch PetscCall(ISDestroy(&perm)); 739566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat)); 749566063dSJacob Faibussowitsch PetscCall(MatDestroy(&fact)); 759566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 769566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 77b122ec5aSJacob Faibussowitsch return 0; 78c4762a1bSJed Brown } 79