xref: /petsc/src/mat/tests/ex117.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
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 
9*9371c9d4SSatish Balay int main(int argc, char **args) {
10c4762a1bSJed Brown   Mat           mat, fact, B;
11c4762a1bSJed Brown   PetscInt      ind1[2], ind2[2];
12c4762a1bSJed Brown   PetscScalar   temp[4];
13c4762a1bSJed Brown   PetscInt      nnz[3];
14c4762a1bSJed Brown   IS            perm, colp;
15c4762a1bSJed Brown   MatFactorInfo info;
16c4762a1bSJed Brown   PetscMPIInt   size;
17c4762a1bSJed Brown 
18327415f7SBarry Smith   PetscFunctionBeginUser;
199566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, 0, help));
209566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
21be096a46SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
22c4762a1bSJed Brown 
23*9371c9d4SSatish Balay   nnz[0] = 2;
24*9371c9d4SSatish Balay   nnz[1] = 1;
25*9371c9d4SSatish Balay   nnz[2] = 1;
269566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqSBAIJ(PETSC_COMM_SELF, 2, 6, 6, 0, nnz, &mat));
27c4762a1bSJed Brown 
28*9371c9d4SSatish Balay   ind1[0] = 0;
29*9371c9d4SSatish Balay   ind1[1] = 1;
30*9371c9d4SSatish Balay   temp[0] = 3;
31*9371c9d4SSatish Balay   temp[1] = 2;
32*9371c9d4SSatish Balay   temp[2] = 0;
33*9371c9d4SSatish Balay   temp[3] = 3;
349566063dSJacob Faibussowitsch   PetscCall(MatSetValues(mat, 2, ind1, 2, ind1, temp, INSERT_VALUES));
35*9371c9d4SSatish Balay   ind2[0] = 4;
36*9371c9d4SSatish Balay   ind2[1] = 5;
37*9371c9d4SSatish Balay   temp[0] = 1;
38*9371c9d4SSatish Balay   temp[1] = 1;
39*9371c9d4SSatish Balay   temp[2] = 2;
40*9371c9d4SSatish Balay   temp[3] = 1;
419566063dSJacob Faibussowitsch   PetscCall(MatSetValues(mat, 2, ind1, 2, ind2, temp, INSERT_VALUES));
42*9371c9d4SSatish Balay   ind1[0] = 2;
43*9371c9d4SSatish Balay   ind1[1] = 3;
44*9371c9d4SSatish Balay   temp[0] = 4;
45*9371c9d4SSatish Balay   temp[1] = 1;
46*9371c9d4SSatish Balay   temp[2] = 1;
47*9371c9d4SSatish Balay   temp[3] = 5;
489566063dSJacob Faibussowitsch   PetscCall(MatSetValues(mat, 2, ind1, 2, ind1, temp, INSERT_VALUES));
49*9371c9d4SSatish Balay   ind1[0] = 4;
50*9371c9d4SSatish Balay   ind1[1] = 5;
51*9371c9d4SSatish Balay   temp[0] = 5;
52*9371c9d4SSatish Balay   temp[1] = 1;
53*9371c9d4SSatish Balay   temp[2] = 1;
54*9371c9d4SSatish Balay   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(MatDuplicate(mat, MAT_SHARE_NONZERO_PATTERN, &B));
61*9371c9d4SSatish Balay   ind1[0] = 0;
62*9371c9d4SSatish Balay   ind1[1] = 1;
63*9371c9d4SSatish Balay   temp[0] = 3;
64*9371c9d4SSatish Balay   temp[1] = 2;
65*9371c9d4SSatish Balay   temp[2] = 0;
66*9371c9d4SSatish Balay   temp[3] = 3;
679566063dSJacob Faibussowitsch   PetscCall(MatSetValues(mat, 2, ind1, 2, ind1, temp, INSERT_VALUES));
68*9371c9d4SSatish Balay   ind2[0] = 4;
69*9371c9d4SSatish Balay   ind2[1] = 5;
70*9371c9d4SSatish Balay   temp[0] = 1;
71*9371c9d4SSatish Balay   temp[1] = 1;
72*9371c9d4SSatish Balay   temp[2] = 2;
73*9371c9d4SSatish Balay   temp[3] = 1;
749566063dSJacob Faibussowitsch   PetscCall(MatSetValues(mat, 2, ind1, 2, ind2, temp, INSERT_VALUES));
75*9371c9d4SSatish Balay   ind1[0] = 2;
76*9371c9d4SSatish Balay   ind1[1] = 3;
77*9371c9d4SSatish Balay   temp[0] = 4;
78*9371c9d4SSatish Balay   temp[1] = 1;
79*9371c9d4SSatish Balay   temp[2] = 1;
80*9371c9d4SSatish Balay   temp[3] = 5;
819566063dSJacob Faibussowitsch   PetscCall(MatSetValues(mat, 2, ind1, 2, ind1, temp, INSERT_VALUES));
82*9371c9d4SSatish Balay   ind1[0] = 4;
83*9371c9d4SSatish Balay   ind1[1] = 5;
84*9371c9d4SSatish Balay   temp[0] = 5;
85*9371c9d4SSatish Balay   temp[1] = 1;
86*9371c9d4SSatish Balay   temp[2] = 1;
87*9371c9d4SSatish Balay   temp[3] = 6;
889566063dSJacob Faibussowitsch   PetscCall(MatSetValues(mat, 2, ind1, 2, ind1, temp, INSERT_VALUES));
89c4762a1bSJed Brown 
909566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY));
919566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY));
92c4762a1bSJed Brown 
939566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "mat: \n"));
949566063dSJacob Faibussowitsch   PetscCall(MatView(mat, PETSC_VIEWER_STDOUT_SELF));
95c4762a1bSJed Brown 
96c4762a1bSJed Brown   /* begin cholesky factorization */
979566063dSJacob Faibussowitsch   PetscCall(MatGetOrdering(mat, MATORDERINGNATURAL, &perm, &colp));
989566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&colp));
99c4762a1bSJed Brown 
100c4762a1bSJed Brown   info.fill = 1.0;
1019566063dSJacob Faibussowitsch   PetscCall(MatGetFactor(mat, MATSOLVERPETSC, MAT_FACTOR_CHOLESKY, &fact));
1029566063dSJacob Faibussowitsch   PetscCall(MatCholeskyFactorSymbolic(fact, mat, perm, &info));
1039566063dSJacob Faibussowitsch   PetscCall(MatCholeskyFactorNumeric(fact, mat, &info));
104c4762a1bSJed Brown 
1059566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&perm));
1069566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mat));
1079566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&fact));
1089566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
1099566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
110b122ec5aSJacob Faibussowitsch   return 0;
111c4762a1bSJed Brown }
112