xref: /petsc/src/mat/tests/ex117.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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   PetscErrorCode ierr;
12c4762a1bSJed Brown   Mat            mat,fact,B;
13c4762a1bSJed Brown   PetscInt       ind1[2],ind2[2];
14c4762a1bSJed Brown   PetscScalar    temp[4];
15c4762a1bSJed Brown   PetscInt       nnz[3];
16c4762a1bSJed Brown   IS             perm,colp;
17c4762a1bSJed Brown   MatFactorInfo  info;
18c4762a1bSJed Brown   PetscMPIInt    size;
19c4762a1bSJed Brown 
20c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,0,help);if (ierr) return ierr;
21*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
222c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!");
23c4762a1bSJed Brown 
24c4762a1bSJed Brown   nnz[0]=2;nnz[1]=1;nnz[2]=1;
25*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
29*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
32*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
35*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
38*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(mat,2,ind1,2,ind1,temp,INSERT_VALUES));
39c4762a1bSJed Brown 
40*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY));
41*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY));
42c4762a1bSJed Brown 
43*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
46*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
49*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
52*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
55*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(mat,2,ind1,2,ind1,temp,INSERT_VALUES));
56c4762a1bSJed Brown 
57*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY));
58*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY));
59c4762a1bSJed Brown 
60*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"mat: \n"));
61*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(mat,PETSC_VIEWER_STDOUT_SELF));
62c4762a1bSJed Brown 
63c4762a1bSJed Brown   /* begin cholesky factorization */
64*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOrdering(mat,MATORDERINGNATURAL,&perm,&colp));
65*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&colp));
66c4762a1bSJed Brown 
67c4762a1bSJed Brown   info.fill=1.0;
68*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetFactor(mat,MATSOLVERPETSC,MAT_FACTOR_CHOLESKY,&fact));
69*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCholeskyFactorSymbolic(fact,mat,perm,&info));
70*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCholeskyFactorNumeric(fact,mat,&info));
71c4762a1bSJed Brown 
72*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&perm));
73*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&mat));
74*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&fact));
75*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
76c4762a1bSJed Brown   ierr = PetscFinalize();
77c4762a1bSJed Brown   return ierr;
78c4762a1bSJed Brown }
79