xref: /petsc/src/mat/tests/ex117.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,0,help));
205f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
212c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!");
22c4762a1bSJed Brown 
23c4762a1bSJed Brown   nnz[0]=2;nnz[1]=1;nnz[2]=1;
245f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSeqSBAIJ(PETSC_COMM_SELF,2,6,6,0,nnz,&mat));
25c4762a1bSJed Brown 
26c4762a1bSJed Brown   ind1[0]=0;ind1[1]=1;
27c4762a1bSJed Brown   temp[0]=3;temp[1]=2;temp[2]=0;temp[3]=3;
285f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(mat,2,ind1,2,ind1,temp,INSERT_VALUES));
29c4762a1bSJed Brown   ind2[0]=4;ind2[1]=5;
30c4762a1bSJed Brown   temp[0]=1;temp[1]=1;temp[2]=2;temp[3]=1;
315f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(mat,2,ind1,2,ind2,temp,INSERT_VALUES));
32c4762a1bSJed Brown   ind1[0]=2;ind1[1]=3;
33c4762a1bSJed Brown   temp[0]=4;temp[1]=1;temp[2]=1;temp[3]=5;
345f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(mat,2,ind1,2,ind1,temp,INSERT_VALUES));
35c4762a1bSJed Brown   ind1[0]=4;ind1[1]=5;
36c4762a1bSJed Brown   temp[0]=5;temp[1]=1;temp[2]=1;temp[3]=6;
375f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(mat,2,ind1,2,ind1,temp,INSERT_VALUES));
38c4762a1bSJed Brown 
395f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY));
405f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY));
41c4762a1bSJed Brown 
425f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDuplicate(mat,MAT_SHARE_NONZERO_PATTERN,&B));
43c4762a1bSJed Brown   ind1[0]=0;ind1[1]=1;
44c4762a1bSJed Brown   temp[0]=3;temp[1]=2;temp[2]=0;temp[3]=3;
455f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(mat,2,ind1,2,ind1,temp,INSERT_VALUES));
46c4762a1bSJed Brown   ind2[0]=4;ind2[1]=5;
47c4762a1bSJed Brown   temp[0]=1;temp[1]=1;temp[2]=2;temp[3]=1;
485f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(mat,2,ind1,2,ind2,temp,INSERT_VALUES));
49c4762a1bSJed Brown   ind1[0]=2;ind1[1]=3;
50c4762a1bSJed Brown   temp[0]=4;temp[1]=1;temp[2]=1;temp[3]=5;
515f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(mat,2,ind1,2,ind1,temp,INSERT_VALUES));
52c4762a1bSJed Brown   ind1[0]=4;ind1[1]=5;
53c4762a1bSJed Brown   temp[0]=5;temp[1]=1;temp[2]=1;temp[3]=6;
545f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(mat,2,ind1,2,ind1,temp,INSERT_VALUES));
55c4762a1bSJed Brown 
565f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY));
575f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY));
58c4762a1bSJed Brown 
595f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"mat: \n"));
605f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(mat,PETSC_VIEWER_STDOUT_SELF));
61c4762a1bSJed Brown 
62c4762a1bSJed Brown   /* begin cholesky factorization */
635f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOrdering(mat,MATORDERINGNATURAL,&perm,&colp));
645f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&colp));
65c4762a1bSJed Brown 
66c4762a1bSJed Brown   info.fill=1.0;
675f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetFactor(mat,MATSOLVERPETSC,MAT_FACTOR_CHOLESKY,&fact));
685f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCholeskyFactorSymbolic(fact,mat,perm,&info));
695f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCholeskyFactorNumeric(fact,mat,&info));
70c4762a1bSJed Brown 
715f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&perm));
725f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&mat));
735f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&fact));
745f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
75*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
76*b122ec5aSJacob Faibussowitsch   return 0;
77c4762a1bSJed Brown }
78